ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_charge_charge.F90
Revision: 1160
Committed: Tue May 11 21:31:15 2004 UTC (20 years, 4 months ago) by gezelter
File size: 4068 byte(s)
Log Message:
Fortran-side changes for group-based cutoffs

File Contents

# User Rev Content
1 gezelter 940 module charge_charge
2    
3     use force_globals
4     use definitions
5     use atype_module
6     use vector_class
7     use simulation
8     use status
9     #ifdef IS_MPI
10     use mpiSimulation
11     #endif
12     implicit none
13    
14     PRIVATE
15 gezelter 1160
16     real(kind=dp), parameter :: pre = 332.06508_DP
17 gezelter 940 logical, save :: haveChargeMap = .false.
18    
19     public::do_charge_pair
20    
21     type :: ChargeList
22     real(kind=DP) :: charge = 0.0_DP
23     end type ChargeList
24    
25     type(ChargeList), dimension(:), allocatable :: ChargeMap
26    
27     contains
28 gezelter 1160
29 gezelter 940 subroutine createChargeMap(status)
30     integer :: nAtypes
31     integer :: status
32     integer :: i
33     real (kind=DP) :: thisCharge
34     logical :: thisProperty
35    
36     status = 0
37    
38     nAtypes = getSize(atypes)
39    
40     if (nAtypes == 0) then
41     status = -1
42     return
43     end if
44    
45     if (.not. allocated(ChargeMap)) then
46     allocate(ChargeMap(nAtypes))
47     endif
48    
49     do i = 1, nAtypes
50    
51     call getElementProperty(atypes, i, "is_Charge", thisProperty)
52    
53     if (thisProperty) then
54     call getElementProperty(atypes, i, "charge", thisCharge)
55     ChargeMap(i)%charge = thisCharge
56     endif
57    
58     end do
59    
60     haveChargeMap = .true.
61    
62     end subroutine createChargeMap
63 gezelter 1160
64 gezelter 1150 subroutine do_charge_pair(atom1, atom2, d, rij, r2, sw, vpair, &
65     pot, f, do_pot, do_stress)
66 gezelter 940
67 gezelter 1150 logical :: do_pot, do_stress
68 gezelter 940
69     integer atom1, atom2, me1, me2, id1, id2
70     integer :: localError
71 gezelter 1150 real(kind=dp) :: rij, r2, q1, q2, sw, vpair
72 gezelter 940 real(kind=dp) :: drdx, drdy, drdz, dudr, fx, fy, fz
73 gezelter 1150 real(kind=dp) :: vterm
74 gezelter 940
75     real( kind = dp ) :: pot
76 gezelter 1150 real( kind = dp ), dimension(3) :: d
77 gezelter 940 real( kind = dp ), dimension(3,nLocal) :: f
78    
79     if (.not.haveChargeMap) then
80     localError = 0
81     call createChargeMap(localError)
82     if ( localError .ne. 0 ) then
83     call handleError("charge-charge", "ChargeMap creation failed!")
84     return
85     end if
86     endif
87    
88     #ifdef IS_MPI
89     me1 = atid_Row(atom1)
90     me2 = atid_Col(atom2)
91     #else
92     me1 = atid(atom1)
93     me2 = atid(atom2)
94     #endif
95    
96     q1 = ChargeMap(me1)%charge
97     q2 = ChargeMap(me2)%charge
98 gezelter 1138
99 gezelter 1160 vterm = pre * q1 * q2 / rij
100 gezelter 1138
101 gezelter 1160 dudr = -sw * vterm / rij
102    
103 gezelter 1150 drdx = d(1) / rij
104     drdy = d(2) / rij
105     drdz = d(3) / rij
106    
107     fx = dudr * drdx
108     fy = dudr * drdy
109     fz = dudr * drdz
110 gezelter 940
111 gezelter 1150 vpair = vpair + vterm
112    
113 gezelter 940 #ifdef IS_MPI
114 gezelter 1150 if (do_pot) then
115 gezelter 1160 pot_Row(atom1) = pot_Row(atom1) + sw*vterm*0.5
116     pot_Col(atom2) = pot_Col(atom2) + sw*vterm*0.5
117 gezelter 1150 endif
118    
119     f_Row(1,atom1) = f_Row(1,atom1) + fx
120     f_Row(2,atom1) = f_Row(2,atom1) + fy
121     f_Row(3,atom1) = f_Row(3,atom1) + fz
122    
123     f_Col(1,atom2) = f_Col(1,atom2) - fx
124     f_Col(2,atom2) = f_Col(2,atom2) - fy
125     f_Col(3,atom2) = f_Col(3,atom2) - fz
126    
127 gezelter 940 #else
128    
129 gezelter 1160 if (do_pot) pot = pot + sw*vterm
130 gezelter 1150
131     f(1,atom1) = f(1,atom1) + fx
132     f(2,atom1) = f(2,atom1) + fy
133     f(3,atom1) = f(3,atom1) + fz
134    
135     f(1,atom2) = f(1,atom2) - fx
136     f(2,atom2) = f(2,atom2) - fy
137     f(3,atom2) = f(3,atom2) - fz
138    
139 gezelter 940 #endif
140 gezelter 1150
141     if (do_stress) then
142    
143 gezelter 940 #ifdef IS_MPI
144 gezelter 1150 id1 = tagRow(atom1)
145     id2 = tagColumn(atom2)
146 gezelter 940 #else
147 gezelter 1150 id1 = atom1
148     id2 = atom2
149 gezelter 940 #endif
150 gezelter 1150
151    
152     if (molMembershipList(id1) .ne. molMembershipList(id2)) then
153    
154     ! because the d vector is the rj - ri vector, and
155     ! because fx, fy, fz are the force on atom i, we need a
156     ! negative sign here:
157    
158     tau_Temp(1) = tau_Temp(1) - d(1) * fx
159     tau_Temp(2) = tau_Temp(2) - d(1) * fy
160     tau_Temp(3) = tau_Temp(3) - d(1) * fz
161     tau_Temp(4) = tau_Temp(4) - d(2) * fx
162     tau_Temp(5) = tau_Temp(5) - d(2) * fy
163     tau_Temp(6) = tau_Temp(6) - d(2) * fz
164     tau_Temp(8) = tau_Temp(8) - d(3) * fy
165     tau_Temp(9) = tau_Temp(9) - d(3) * fz
166    
167    
168     virial_Temp = virial_Temp + &
169     (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
170    
171 gezelter 940 endif
172     endif
173 gezelter 1150
174    
175 gezelter 940 return
176     end subroutine do_charge_pair
177    
178     end module charge_charge