1 |
gezelter |
1930 |
!! |
2 |
|
|
!! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. |
3 |
|
|
!! |
4 |
|
|
!! The University of Notre Dame grants you ("Licensee") a |
5 |
|
|
!! non-exclusive, royalty free, license to use, modify and |
6 |
|
|
!! redistribute this software in source and binary code form, provided |
7 |
|
|
!! that the following conditions are met: |
8 |
|
|
!! |
9 |
|
|
!! 1. Acknowledgement of the program authors must be made in any |
10 |
|
|
!! publication of scientific results based in part on use of the |
11 |
|
|
!! program. An acceptable form of acknowledgement is citation of |
12 |
|
|
!! the article in which the program was described (Matthew |
13 |
|
|
!! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher |
14 |
|
|
!! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented |
15 |
|
|
!! Parallel Simulation Engine for Molecular Dynamics," |
16 |
|
|
!! J. Comput. Chem. 26, pp. 252-271 (2005)) |
17 |
|
|
!! |
18 |
|
|
!! 2. Redistributions of source code must retain the above copyright |
19 |
|
|
!! notice, this list of conditions and the following disclaimer. |
20 |
|
|
!! |
21 |
|
|
!! 3. Redistributions in binary form must reproduce the above copyright |
22 |
|
|
!! notice, this list of conditions and the following disclaimer in the |
23 |
|
|
!! documentation and/or other materials provided with the |
24 |
|
|
!! distribution. |
25 |
|
|
!! |
26 |
|
|
!! This software is provided "AS IS," without a warranty of any |
27 |
|
|
!! kind. All express or implied conditions, representations and |
28 |
|
|
!! warranties, including any implied warranty of merchantability, |
29 |
|
|
!! fitness for a particular purpose or non-infringement, are hereby |
30 |
|
|
!! excluded. The University of Notre Dame and its licensors shall not |
31 |
|
|
!! be liable for any damages suffered by licensee as a result of |
32 |
|
|
!! using, modifying or distributing the software or its |
33 |
|
|
!! derivatives. In no event will the University of Notre Dame or its |
34 |
|
|
!! licensors be liable for any lost revenue, profit or data, or for |
35 |
|
|
!! direct, indirect, special, consequential, incidental or punitive |
36 |
|
|
!! damages, however caused and regardless of the theory of liability, |
37 |
|
|
!! arising out of the use of or inability to use software, even if the |
38 |
|
|
!! University of Notre Dame has been advised of the possibility of |
39 |
|
|
!! such damages. |
40 |
|
|
!! |
41 |
|
|
|
42 |
gezelter |
1608 |
module reaction_field |
43 |
|
|
use force_globals |
44 |
|
|
use definitions |
45 |
|
|
use atype_module |
46 |
|
|
use vector_class |
47 |
|
|
use simulation |
48 |
|
|
use status |
49 |
|
|
#ifdef IS_MPI |
50 |
|
|
use mpiSimulation |
51 |
|
|
#endif |
52 |
|
|
implicit none |
53 |
|
|
|
54 |
|
|
PRIVATE |
55 |
|
|
|
56 |
|
|
real(kind=dp), save :: rrf = 1.0_dp |
57 |
|
|
real(kind=dp), save :: rt |
58 |
|
|
real(kind=dp), save :: dielect = 1.0_dp |
59 |
|
|
real(kind=dp), save :: rrfsq = 1.0_dp |
60 |
|
|
real(kind=dp), save :: pre |
61 |
|
|
logical, save :: haveCutoffs = .false. |
62 |
|
|
logical, save :: haveMomentMap = .false. |
63 |
|
|
logical, save :: haveDielectric = .false. |
64 |
|
|
|
65 |
|
|
type :: MomentList |
66 |
|
|
real(kind=DP) :: dipole_moment = 0.0_DP |
67 |
|
|
end type MomentList |
68 |
|
|
|
69 |
|
|
type(MomentList), dimension(:),allocatable :: MomentMap |
70 |
|
|
|
71 |
|
|
PUBLIC::initialize_rf |
72 |
|
|
PUBLIC::setCutoffsRF |
73 |
|
|
PUBLIC::accumulate_rf |
74 |
|
|
PUBLIC::accumulate_self_rf |
75 |
|
|
PUBLIC::reaction_field_final |
76 |
|
|
PUBLIC::rf_correct_forces |
77 |
|
|
|
78 |
|
|
contains |
79 |
|
|
|
80 |
|
|
subroutine initialize_rf(this_dielect) |
81 |
|
|
real(kind=dp), intent(in) :: this_dielect |
82 |
|
|
|
83 |
|
|
dielect = this_dielect |
84 |
|
|
|
85 |
|
|
pre = 14.38362d0*2.0d0*(dielect-1.0d0)/((2.0d0*dielect+1.0d0)*rrfsq*rrf) |
86 |
|
|
|
87 |
|
|
haveDielectric = .true. |
88 |
|
|
|
89 |
|
|
return |
90 |
|
|
end subroutine initialize_rf |
91 |
|
|
|
92 |
|
|
subroutine setCutoffsRF( this_rrf, this_rt ) |
93 |
|
|
|
94 |
|
|
real(kind=dp), intent(in) :: this_rrf, this_rt |
95 |
|
|
|
96 |
|
|
rrf = this_rrf |
97 |
|
|
rt = this_rt |
98 |
|
|
|
99 |
|
|
rrfsq = rrf * rrf |
100 |
|
|
pre = 14.38362d0*2.0d0*(dielect-1.0d0)/((2.0d0*dielect+1.0d0)*rrfsq*rrf) |
101 |
|
|
|
102 |
|
|
haveCutoffs = .true. |
103 |
|
|
|
104 |
|
|
end subroutine setCutoffsRF |
105 |
|
|
|
106 |
|
|
subroutine createMomentMap(status) |
107 |
|
|
integer :: nAtypes |
108 |
|
|
integer :: status |
109 |
|
|
integer :: i |
110 |
|
|
real (kind=DP) :: thisDP |
111 |
|
|
logical :: thisProperty |
112 |
|
|
|
113 |
|
|
status = 0 |
114 |
|
|
|
115 |
|
|
nAtypes = getSize(atypes) |
116 |
|
|
|
117 |
|
|
if (nAtypes == 0) then |
118 |
|
|
status = -1 |
119 |
|
|
return |
120 |
|
|
end if |
121 |
|
|
|
122 |
|
|
if (.not. allocated(MomentMap)) then |
123 |
|
|
allocate(MomentMap(nAtypes)) |
124 |
|
|
endif |
125 |
|
|
|
126 |
tim |
2058 |
!!do i = 1, nAtypes |
127 |
|
|
!! |
128 |
|
|
!! call getElementProperty(atypes, i, "is_Dipole", thisProperty) |
129 |
|
|
!! |
130 |
|
|
!! if (thisProperty) then |
131 |
|
|
!! call getElementProperty(atypes, i, "dipole_moment", thisDP) |
132 |
|
|
!! MomentMap(i)%dipole_moment = thisDP |
133 |
|
|
!! endif |
134 |
|
|
!! |
135 |
|
|
!!end do |
136 |
gezelter |
1608 |
|
137 |
|
|
haveMomentMap = .true. |
138 |
|
|
|
139 |
|
|
end subroutine createMomentMap |
140 |
|
|
|
141 |
gezelter |
1930 |
subroutine accumulate_rf(atom1, atom2, rij, eFrame, taper) |
142 |
gezelter |
1608 |
|
143 |
|
|
integer, intent(in) :: atom1, atom2 |
144 |
|
|
real (kind = dp), intent(in) :: rij |
145 |
gezelter |
1930 |
real (kind = dp), dimension(9,nLocal) :: eFrame |
146 |
gezelter |
1608 |
|
147 |
|
|
integer :: me1, me2 |
148 |
|
|
real (kind = dp), intent(in) :: taper |
149 |
|
|
real (kind = dp):: mu1, mu2 |
150 |
|
|
real (kind = dp), dimension(3) :: ul1 |
151 |
|
|
real (kind = dp), dimension(3) :: ul2 |
152 |
|
|
|
153 |
|
|
integer :: localError |
154 |
|
|
|
155 |
|
|
if ((.not.haveDielectric).or.(.not.haveCutoffs)) then |
156 |
|
|
write(default_error,*) 'Reaction field not initialized!' |
157 |
|
|
return |
158 |
|
|
endif |
159 |
|
|
|
160 |
|
|
if (.not.haveMomentMap) then |
161 |
|
|
localError = 0 |
162 |
|
|
call createMomentMap(localError) |
163 |
|
|
if ( localError .ne. 0 ) then |
164 |
|
|
call handleError("reaction-field", "MomentMap creation failed!") |
165 |
|
|
return |
166 |
|
|
end if |
167 |
|
|
endif |
168 |
|
|
|
169 |
|
|
|
170 |
|
|
#ifdef IS_MPI |
171 |
|
|
me1 = atid_Row(atom1) |
172 |
gezelter |
1930 |
ul1(1) = eFrame_Row(3,atom1) |
173 |
|
|
ul1(2) = eFrame_Row(6,atom1) |
174 |
|
|
ul1(3) = eFrame_Row(9,atom1) |
175 |
gezelter |
1608 |
|
176 |
|
|
me2 = atid_Col(atom2) |
177 |
gezelter |
1930 |
ul2(1) = eFrame_Col(3,atom2) |
178 |
|
|
ul2(2) = eFrame_Col(6,atom2) |
179 |
|
|
ul2(3) = eFrame_Col(9,atom2) |
180 |
gezelter |
1608 |
#else |
181 |
|
|
me1 = atid(atom1) |
182 |
gezelter |
1930 |
ul1(1) = eFrame(3,atom1) |
183 |
|
|
ul1(2) = eFrame(6,atom1) |
184 |
|
|
ul1(3) = eFrame(9,atom1) |
185 |
gezelter |
1608 |
|
186 |
|
|
me2 = atid(atom2) |
187 |
gezelter |
1930 |
ul2(1) = eFrame(3,atom2) |
188 |
|
|
ul2(2) = eFrame(6,atom2) |
189 |
|
|
ul2(3) = eFrame(9,atom2) |
190 |
gezelter |
1608 |
#endif |
191 |
|
|
|
192 |
|
|
mu1 = MomentMap(me1)%dipole_moment |
193 |
|
|
mu2 = MomentMap(me2)%dipole_moment |
194 |
|
|
|
195 |
|
|
#ifdef IS_MPI |
196 |
|
|
rf_Row(1,atom1) = rf_Row(1,atom1) + ul2(1)*mu2*taper |
197 |
|
|
rf_Row(2,atom1) = rf_Row(2,atom1) + ul2(2)*mu2*taper |
198 |
|
|
rf_Row(3,atom1) = rf_Row(3,atom1) + ul2(3)*mu2*taper |
199 |
|
|
|
200 |
|
|
rf_Col(1,atom2) = rf_Col(1,atom2) + ul1(1)*mu1*taper |
201 |
|
|
rf_Col(2,atom2) = rf_Col(2,atom2) + ul1(2)*mu1*taper |
202 |
|
|
rf_Col(3,atom2) = rf_Col(3,atom2) + ul1(3)*mu1*taper |
203 |
|
|
#else |
204 |
|
|
rf(1,atom1) = rf(1,atom1) + ul2(1)*mu2*taper |
205 |
|
|
rf(2,atom1) = rf(2,atom1) + ul2(2)*mu2*taper |
206 |
|
|
rf(3,atom1) = rf(3,atom1) + ul2(3)*mu2*taper |
207 |
|
|
|
208 |
|
|
rf(1,atom2) = rf(1,atom2) + ul1(1)*mu1*taper |
209 |
|
|
rf(2,atom2) = rf(2,atom2) + ul1(2)*mu1*taper |
210 |
|
|
rf(3,atom2) = rf(3,atom2) + ul1(3)*mu1*taper |
211 |
|
|
#endif |
212 |
|
|
|
213 |
|
|
|
214 |
|
|
return |
215 |
|
|
end subroutine accumulate_rf |
216 |
|
|
|
217 |
gezelter |
1930 |
subroutine accumulate_self_rf(atom1, mu1, eFrame) |
218 |
gezelter |
1608 |
|
219 |
|
|
integer, intent(in) :: atom1 |
220 |
|
|
real(kind=dp), intent(in) :: mu1 |
221 |
gezelter |
1930 |
real(kind=dp), dimension(9,nLocal) :: eFrame |
222 |
gezelter |
1608 |
|
223 |
|
|
!! should work for both MPI and non-MPI version since this is not pairwise. |
224 |
gezelter |
1930 |
rf(1,atom1) = rf(1,atom1) + eFrame(3,atom1)*mu1 |
225 |
|
|
rf(2,atom1) = rf(2,atom1) + eFrame(6,atom1)*mu1 |
226 |
|
|
rf(3,atom1) = rf(3,atom1) + eFrame(9,atom1)*mu1 |
227 |
gezelter |
1608 |
|
228 |
|
|
return |
229 |
|
|
end subroutine accumulate_self_rf |
230 |
|
|
|
231 |
gezelter |
1930 |
subroutine reaction_field_final(a1, mu1, eFrame, rfpot, t, do_pot) |
232 |
gezelter |
1608 |
|
233 |
|
|
integer, intent(in) :: a1 |
234 |
|
|
real (kind=dp), intent(in) :: mu1 |
235 |
|
|
real (kind=dp), intent(inout) :: rfpot |
236 |
|
|
logical, intent(in) :: do_pot |
237 |
gezelter |
1930 |
real (kind = dp), dimension(9,nLocal) :: eFrame |
238 |
gezelter |
1608 |
real (kind = dp), dimension(3,nLocal) :: t |
239 |
|
|
|
240 |
|
|
integer :: localError |
241 |
|
|
|
242 |
|
|
if ((.not.haveDielectric).or.(.not.haveCutoffs)) then |
243 |
|
|
write(default_error,*) 'Reaction field not initialized!' |
244 |
|
|
return |
245 |
|
|
endif |
246 |
|
|
|
247 |
|
|
if (.not.haveMomentMap) then |
248 |
|
|
localError = 0 |
249 |
|
|
call createMomentMap(localError) |
250 |
|
|
if ( localError .ne. 0 ) then |
251 |
|
|
call handleError("reaction-field", "MomentMap creation failed!") |
252 |
|
|
return |
253 |
|
|
end if |
254 |
|
|
endif |
255 |
|
|
|
256 |
|
|
! compute torques on dipoles: |
257 |
|
|
! pre converts from mu in units of debye to kcal/mol |
258 |
|
|
|
259 |
|
|
! The torque contribution is dipole cross reaction_field |
260 |
|
|
|
261 |
gezelter |
1930 |
t(1,a1) = t(1,a1) + pre*mu1*(eFrame(6,a1)*rf(3,a1) - eFrame(9,a1)*rf(2,a1)) |
262 |
|
|
t(2,a1) = t(2,a1) + pre*mu1*(eFrame(9,a1)*rf(1,a1) - eFrame(3,a1)*rf(3,a1)) |
263 |
|
|
t(3,a1) = t(3,a1) + pre*mu1*(eFrame(3,a1)*rf(2,a1) - eFrame(6,a1)*rf(1,a1)) |
264 |
gezelter |
1608 |
|
265 |
|
|
! the potential contribution is -1/2 dipole dot reaction_field |
266 |
|
|
|
267 |
|
|
if (do_pot) then |
268 |
|
|
rfpot = rfpot - 0.5d0 * pre * mu1 * & |
269 |
gezelter |
1930 |
(rf(1,a1)*eFrame(3,a1) + rf(2,a1)*eFrame(6,a1) + rf(3,a1)*eFrame(9,a1)) |
270 |
gezelter |
1608 |
endif |
271 |
|
|
|
272 |
|
|
return |
273 |
|
|
end subroutine reaction_field_final |
274 |
|
|
|
275 |
gezelter |
1930 |
subroutine rf_correct_forces(atom1, atom2, d, rij, eFrame, taper, f, fpair) |
276 |
gezelter |
1608 |
|
277 |
|
|
integer, intent(in) :: atom1, atom2 |
278 |
|
|
real(kind=dp), dimension(3), intent(in) :: d |
279 |
|
|
real(kind=dp), intent(in) :: rij, taper |
280 |
gezelter |
1930 |
real( kind = dp ), dimension(9,nLocal) :: eFrame |
281 |
gezelter |
1608 |
real( kind = dp ), dimension(3,nLocal) :: f |
282 |
|
|
real( kind = dp ), dimension(3), intent(inout) :: fpair |
283 |
|
|
|
284 |
|
|
real (kind = dp), dimension(3) :: ul1 |
285 |
|
|
real (kind = dp), dimension(3) :: ul2 |
286 |
|
|
real (kind = dp) :: dtdr |
287 |
|
|
real (kind = dp) :: dudx, dudy, dudz, u1dotu2 |
288 |
|
|
integer :: me1, me2, id1, id2 |
289 |
|
|
real (kind = dp) :: mu1, mu2 |
290 |
|
|
|
291 |
|
|
integer :: localError |
292 |
|
|
|
293 |
|
|
if ((.not.haveDielectric).or.(.not.haveCutoffs)) then |
294 |
|
|
write(default_error,*) 'Reaction field not initialized!' |
295 |
|
|
return |
296 |
|
|
endif |
297 |
|
|
|
298 |
|
|
if (.not.haveMomentMap) then |
299 |
|
|
localError = 0 |
300 |
|
|
call createMomentMap(localError) |
301 |
|
|
if ( localError .ne. 0 ) then |
302 |
|
|
call handleError("reaction-field", "MomentMap creation failed!") |
303 |
|
|
return |
304 |
|
|
end if |
305 |
|
|
endif |
306 |
|
|
|
307 |
|
|
if (rij.le.rrf) then |
308 |
|
|
|
309 |
|
|
if (rij.lt.rt) then |
310 |
|
|
dtdr = 0.0d0 |
311 |
|
|
else |
312 |
|
|
! write(*,*) 'rf correct in taper region' |
313 |
|
|
dtdr = 6.0d0*(rij*rij - rij*rt - rij*rrf +rrf*rt)/((rrf-rt)**3) |
314 |
|
|
endif |
315 |
|
|
|
316 |
|
|
#ifdef IS_MPI |
317 |
|
|
me1 = atid_Row(atom1) |
318 |
gezelter |
1930 |
ul1(1) = eFrame_Row(3,atom1) |
319 |
|
|
ul1(2) = eFrame_Row(6,atom1) |
320 |
|
|
ul1(3) = eFrame_Row(9,atom1) |
321 |
gezelter |
1608 |
|
322 |
|
|
me2 = atid_Col(atom2) |
323 |
gezelter |
1930 |
ul2(1) = eFrame_Col(3,atom2) |
324 |
|
|
ul2(2) = eFrame_Col(6,atom2) |
325 |
|
|
ul2(3) = eFrame_Col(9,atom2) |
326 |
gezelter |
1608 |
#else |
327 |
|
|
me1 = atid(atom1) |
328 |
gezelter |
1930 |
ul1(1) = eFrame(3,atom1) |
329 |
|
|
ul1(2) = eFrame(6,atom1) |
330 |
|
|
ul1(3) = eFrame(9,atom1) |
331 |
gezelter |
1608 |
|
332 |
|
|
me2 = atid(atom2) |
333 |
gezelter |
1930 |
ul2(1) = eFrame(3,atom2) |
334 |
|
|
ul2(2) = eFrame(6,atom2) |
335 |
|
|
ul2(3) = eFrame(9,atom2) |
336 |
gezelter |
1608 |
#endif |
337 |
|
|
|
338 |
|
|
mu1 = MomentMap(me1)%dipole_moment |
339 |
|
|
mu2 = MomentMap(me2)%dipole_moment |
340 |
|
|
|
341 |
|
|
u1dotu2 = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3) |
342 |
|
|
|
343 |
|
|
dudx = - pre*mu1*mu2*u1dotu2*dtdr*d(1)/rij |
344 |
|
|
dudy = - pre*mu1*mu2*u1dotu2*dtdr*d(2)/rij |
345 |
|
|
dudz = - pre*mu1*mu2*u1dotu2*dtdr*d(3)/rij |
346 |
|
|
|
347 |
|
|
#ifdef IS_MPI |
348 |
|
|
f_Row(1,atom1) = f_Row(1,atom1) + dudx |
349 |
|
|
f_Row(2,atom1) = f_Row(2,atom1) + dudy |
350 |
|
|
f_Row(3,atom1) = f_Row(3,atom1) + dudz |
351 |
|
|
|
352 |
|
|
f_Col(1,atom2) = f_Col(1,atom2) - dudx |
353 |
|
|
f_Col(2,atom2) = f_Col(2,atom2) - dudy |
354 |
|
|
f_Col(3,atom2) = f_Col(3,atom2) - dudz |
355 |
|
|
#else |
356 |
|
|
f(1,atom1) = f(1,atom1) + dudx |
357 |
|
|
f(2,atom1) = f(2,atom1) + dudy |
358 |
|
|
f(3,atom1) = f(3,atom1) + dudz |
359 |
|
|
|
360 |
|
|
f(1,atom2) = f(1,atom2) - dudx |
361 |
|
|
f(2,atom2) = f(2,atom2) - dudy |
362 |
|
|
f(3,atom2) = f(3,atom2) - dudz |
363 |
|
|
#endif |
364 |
|
|
|
365 |
|
|
#ifdef IS_MPI |
366 |
|
|
id1 = AtomRowToGlobal(atom1) |
367 |
|
|
id2 = AtomColToGlobal(atom2) |
368 |
|
|
#else |
369 |
|
|
id1 = atom1 |
370 |
|
|
id2 = atom2 |
371 |
|
|
#endif |
372 |
|
|
|
373 |
|
|
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
374 |
|
|
|
375 |
|
|
fpair(1) = fpair(1) + dudx |
376 |
|
|
fpair(2) = fpair(2) + dudy |
377 |
|
|
fpair(3) = fpair(3) + dudz |
378 |
|
|
|
379 |
|
|
endif |
380 |
|
|
|
381 |
|
|
end if |
382 |
|
|
return |
383 |
|
|
end subroutine rf_correct_forces |
384 |
|
|
end module reaction_field |