ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/calc_reaction_field.F90
Revision: 330
Committed: Wed Mar 12 23:15:46 2003 UTC (21 years, 4 months ago) by gezelter
File size: 7382 byte(s)
Log Message:
forceGlobals is gone (part3)

File Contents

# Content
1 module reaction_field
2 use simulation
3 use definitions
4 use atype_module
5 use vector_class
6 #ifdef IS_MPI
7 use mpiSimulation
8 #endif
9 implicit none
10
11 real(kind=dp), save :: rrf
12 real(kind=dp), save :: rt
13 real(kind=dp), save :: dielect
14 real(kind=dp), save :: rrfsq
15 real(kind=dp), save :: pre
16 logical, save :: rf_initialized = .false.
17
18 contains
19
20 subroutine initialize_rf()
21 rrf = getRrf()
22 rt = getRt()
23 dielect = getDielect()
24
25 rrfsq = rrf * rrf
26 pre = 14.38362d0*2.0d0*(dielect-1.0d0)/((2.0d0*dielect+1.0d0)*rrfsq*rrf)
27
28 rf_initialized = .true.
29
30 return
31 end subroutine initialize_rf
32
33 subroutine accumulate_rf(atom1, atom2, rij, u_l)
34
35 integer, intent(in) :: atom1, atom2
36 real (kind = dp), intent(in) :: rij
37 real (kind = dp), dimension(3, getNlocal()) :: u_l
38
39 integer :: me1, me2
40 real (kind = dp) :: taper, mu1, mu2
41 real (kind = dp), dimension(3) :: ul1
42 real (kind = dp), dimension(3) :: ul2
43
44 if (.not.rf_initialized) then
45 write(default_error,*) 'Reaction field not initialized!'
46 return
47 endif
48
49 if (rij.le.rrf) then
50
51 if (rij.lt.rt) then
52 taper = 1.0d0
53 else
54 taper = (rrf + 2.0d0*rij - 3.0d0*rt)*(rrf-rij)**2/ ((rrf-rt)**3)
55 endif
56
57 #ifdef IS_MPI
58 me1 = atid_Row(atom1)
59 ul1(1) = u_l_Row(1,atom1)
60 ul1(2) = u_l_Row(2,atom1)
61 ul1(3) = u_l_Row(3,atom1)
62
63 me2 = atid_Col(atom2)
64 ul2(1) = u_l_Col(1,atom2)
65 ul2(2) = u_l_Col(2,atom2)
66 ul2(3) = u_l_Col(3,atom2)
67 #else
68 me1 = atid(atom1)
69 ul1(1) = u_l(1,atom1)
70 ul1(2) = u_l(2,atom1)
71 ul1(3) = u_l(3,atom1)
72
73 me2 = atid(atom2)
74 ul2(1) = u_l(1,atom2)
75 ul2(2) = u_l(2,atom2)
76 ul2(3) = u_l(3,atom2)
77 #endif
78
79 call getElementProperty(atypes, me1, "dipole_moment", mu1)
80 call getElementProperty(atypes, me2, "dipole_moment", mu2)
81
82
83 #ifdef IS_MPI
84 rf_Row(1,atom1) = rf_Row(1,atom1) + ul2(1)*mu2*taper
85 rf_Row(2,atom1) = rf_Row(2,atom1) + ul2(2)*mu2*taper
86 rf_Row(3,atom1) = rf_Row(3,atom1) + ul2(3)*mu2*taper
87
88 rf_Col(1,atom2) = rf_Col(1,atom2) + ul1(1)*mu2*taper
89 rf_Col(2,atom2) = rf_Col(2,atom2) + ul1(2)*mu2*taper
90 rf_Col(3,atom2) = rf_Col(3,atom2) + ul1(3)*mu2*taper
91 #else
92 rf(1,atom1) = rf(1,atom1) + ul2(1)*mu2*taper
93 rf(2,atom1) = rf(2,atom1) + ul2(2)*mu2*taper
94 rf(3,atom1) = rf(3,atom1) + ul2(3)*mu2*taper
95
96 rf(1,atom2) = rf(1,atom2) + ul1(1)*mu2*taper
97 rf(2,atom2) = rf(2,atom2) + ul1(2)*mu2*taper
98 rf(3,atom2) = rf(3,atom2) + ul1(3)*mu2*taper
99 #endif
100
101 endif
102 return
103 end subroutine accumulate_rf
104
105 subroutine accumulate_self_rf(atom1, mu1, u_l)
106
107 integer, intent(in) :: atom1
108 real(kind=dp), intent(in) :: mu1
109 real(kind=dp), dimension(3,getNlocal()) :: u_l
110
111 !! should work for both MPI and non-MPI version since this is not pairwise.
112 rf(1,atom1) = rf(1,atom1) + u_l(1,atom1)*mu1
113 rf(2,atom1) = rf(2,atom1) + u_l(2,atom1)*mu1
114 rf(3,atom1) = rf(3,atom1) + u_l(3,atom1)*mu1
115
116 return
117 end subroutine accumulate_self_rf
118
119 subroutine reaction_field_final(a1, mu1, u_l, rfpot, t, do_pot)
120
121 integer, intent(in) :: a1
122 real (kind=dp), intent(in) :: mu1
123 real (kind=dp), intent(inout) :: rfpot
124 logical, intent(in) :: do_pot
125 real (kind = dp), dimension(3, getNlocal()) :: u_l
126 real (kind = dp), dimension(3, getNlocal()) :: t
127
128 if (.not.rf_initialized) then
129 write(default_error,*) 'Reaction field not initialized!'
130 return
131 endif
132
133 ! compute torques on dipoles:
134 ! pre converts from mu in units of debye to kcal/mol
135
136 ! The torque contribution is dipole cross reaction_field
137
138 t(1,a1) = t(1,a1) + pre*mu1*(u_l(2,a1)*rf(3,a1) - u_l(3,a1)*rf(2,a1))
139 t(2,a1) = t(2,a1) + pre*mu1*(u_l(3,a1)*rf(1,a1) - u_l(1,a1)*rf(3,a1))
140 t(3,a1) = t(3,a1) + pre*mu1*(u_l(1,a1)*rf(2,a1) - u_l(2,a1)*rf(1,a1))
141
142 ! the potential contribution is -1/2 dipole dot reaction_field
143
144 if (do_pot) then
145 rfpot = rfpot - 0.5d0 * pre * mu1 * &
146 (rf(1,a1)*u_l(1,a1) + rf(2,a1)*u_l(2,a1) + rf(3,a1)*u_l(3,a1))
147 endif
148
149 return
150 end subroutine reaction_field_final
151
152 subroutine rf_correct_forces(atom1, atom2, d, rij, u_l, f, do_stress)
153
154 integer, intent(in) :: atom1, atom2
155 real(kind=dp), dimension(3), intent(in) :: d
156 real(kind=dp), intent(in) :: rij
157 real( kind = dp ), dimension(3,getNlocal()) :: u_l
158 real( kind = dp ), dimension(3,getNlocal()) :: f
159 logical, intent(in) :: do_stress
160
161 real (kind = dp), dimension(3) :: ul1
162 real (kind = dp), dimension(3) :: ul2
163 real (kind = dp) :: dtdr
164 real (kind = dp) :: dudx, dudy, dudz, u1dotu2
165 integer :: me1, me2
166 real (kind = dp) :: mu1, mu2
167
168 if (.not.rf_initialized) then
169 write(default_error,*) 'Reaction field not initialized!'
170 return
171 endif
172
173 if (rij.le.rrf) then
174
175 if (rij.lt.rt) then
176 dtdr = 0.0d0
177 else
178 dtdr = 6.0d0*(rij*rij - rij*rt - rij*rrf +rrf*rt)/((rrf-rt)**3)
179 endif
180
181 #ifdef IS_MPI
182 me1 = atid_Row(atom1)
183 ul1(1) = u_l_Row(1,atom1)
184 ul1(2) = u_l_Row(2,atom1)
185 ul1(3) = u_l_Row(3,atom1)
186
187 me2 = atid_Col(atom2)
188 ul2(1) = u_l_Col(1,atom2)
189 ul2(2) = u_l_Col(2,atom2)
190 ul2(3) = u_l_Col(3,atom2)
191 #else
192 me1 = atid(atom1)
193 ul1(1) = u_l(1,atom1)
194 ul1(2) = u_l(2,atom1)
195 ul1(3) = u_l(3,atom1)
196
197 me2 = atid(atom2)
198 ul2(1) = u_l(1,atom2)
199 ul2(2) = u_l(2,atom2)
200 ul2(3) = u_l(3,atom2)
201 #endif
202
203 call getElementProperty(atypes, me1, "dipole_moment", mu1)
204 call getElementProperty(atypes, me2, "dipole_moment", mu2)
205
206 u1dotu2 = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3)
207
208 dudx = - pre*mu1*mu2*u1dotu2*dtdr*d(1)/rij
209 dudy = - pre*mu1*mu2*u1dotu2*dtdr*d(2)/rij
210 dudz = - pre*mu1*mu2*u1dotu2*dtdr*d(3)/rij
211
212 #ifdef IS_MPI
213 f_Row(1,atom1) = f_Row(1,atom1) + dudx
214 f_Row(2,atom1) = f_Row(2,atom1) + dudy
215 f_Row(3,atom1) = f_Row(3,atom1) + dudz
216
217 f_Col(1,atom2) = f_Col(1,atom2) - dudx
218 f_Col(2,atom2) = f_Col(2,atom2) - dudy
219 f_Col(3,atom2) = f_Col(3,atom2) - dudz
220 #else
221 f(1,atom1) = f(1,atom1) + dudx
222 f(2,atom1) = f(2,atom1) + dudy
223 f(3,atom1) = f(3,atom1) + dudz
224
225 f(1,atom2) = f(1,atom2) - dudx
226 f(2,atom2) = f(2,atom2) - dudy
227 f(3,atom2) = f(3,atom2) - dudz
228 #endif
229
230 if (do_stress) then
231 tau_Temp(1) = tau_Temp(1) + dudx * d(1)
232 tau_Temp(2) = tau_Temp(2) + dudx * d(2)
233 tau_Temp(3) = tau_Temp(3) + dudx * d(3)
234 tau_Temp(4) = tau_Temp(4) + dudy * d(1)
235 tau_Temp(5) = tau_Temp(5) + dudy * d(2)
236 tau_Temp(6) = tau_Temp(6) + dudy * d(3)
237 tau_Temp(7) = tau_Temp(7) + dudz * d(1)
238 tau_Temp(8) = tau_Temp(8) + dudz * d(2)
239 tau_Temp(9) = tau_Temp(9) + dudz * d(3)
240 virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
241 endif
242 endif
243
244 return
245 end subroutine rf_correct_forces
246 end module reaction_field