1 |
module long_range_forces |
2 |
use definitions, ONLY : dp, ndim |
3 |
use simulation |
4 |
implicit none |
5 |
|
6 |
|
7 |
|
8 |
|
9 |
real(kind = dp ) :: rlistsq |
10 |
real(kind = dp ) :: rcutsq |
11 |
|
12 |
|
13 |
|
14 |
|
15 |
|
16 |
|
17 |
|
18 |
subroutine long_range_force ( check_exclude, & |
19 |
pair_i, pair_j, n_pairs, natoms, sigma, epslon, & |
20 |
box, & |
21 |
rx, ry, rz, fx, fy, fz, & |
22 |
v, & |
23 |
update, point, list, maxnab, & |
24 |
rx0, ry0, rz0, & |
25 |
ux, uy, uz, & |
26 |
tx, ty, tz, & |
27 |
ex, ey, ez, & |
28 |
mu, rrf, rtaper, & |
29 |
dielectric, & |
30 |
is_dipole, is_vdw, is_lj, is_ssd, & |
31 |
Axx, Axy, Axz, & |
32 |
Ayx, Ayy, Ayz, & |
33 |
Azx, Azy, Azz ) |
34 |
|
35 |
implicit none |
36 |
|
37 |
! This routine implements the Verlet neighbor list for all of |
38 |
! the long range force calculations. After a pair has been |
39 |
! determined to lie within the neighbor list, the pair is passed |
40 |
! to the appropriate force functions. |
41 |
|
42 |
|
43 |
|
44 |
! Passed Arguments |
45 |
|
46 |
integer :: n_pairs ! the number of excluded pairs |
47 |
integer :: natoms ! the number of atoms in the system |
48 |
integer :: maxnab ! the max number of neighbors for the verlet list |
49 |
|
50 |
logical(kind=2) :: check_exclude ! boolean to check for exclusion pairs |
51 |
logical(kind=2) :: update ! boolean to update the neighbor list |
52 |
|
53 |
double precision rcut ! the VDW/LJ cutoff radius |
54 |
double precision rlist ! Verlet list cutoff |
55 |
double precision box_x, box_y, box_z ! periodic boundry conditions |
56 |
double precision rrf ! dipole reaction field cutoff |
57 |
double precision rtaper ! the taper radius of the rxn field |
58 |
double precision dielectric ! the dielectric cnst. of the rxn field |
59 |
double precision v ! the potential energy |
60 |
|
61 |
|
62 |
|
63 |
! Passed Arrays |
64 |
|
65 |
integer, dimension(n_pairs):: pair_i, pair_j ! the excluded pairs i and j |
66 |
integer, dimension(natoms) :: point ! the Verlet list ptr array |
67 |
integer, dimension(maxnab) :: list ! the verlet list |
68 |
|
69 |
logical(kind=2), dimension(natoms) :: is_dipole ! dipole boolean array |
70 |
logical(kind=2), dimension(natoms) :: is_vdw ! VDW boolean array |
71 |
logical(kind=2), dimension(natoms) :: is_lj ! LJ boolean array |
72 |
logical(kind=2), dimension(natoms) :: is_ssd ! ssd boolean array |
73 |
|
74 |
double precision, dimension(natoms) :: sigma ! VDW/LJ distance prmtr. |
75 |
double precision, dimension(natoms) :: epslon ! VDW/LJ well depth prmtr. |
76 |
double precision, dimension(natoms) :: rx, ry, rz ! positions |
77 |
double precision, dimension(natoms) :: fx, fy, fz ! forces |
78 |
double precision, dimension(natoms) :: rx0, ry0, rz0 ! intial verlet positions |
79 |
double precision, dimension(natoms) :: ux, uy, uz ! dipole unit vectors |
80 |
double precision, dimension(natoms) :: tx, ty, tz ! torsions |
81 |
double precision, dimension(natoms) :: ex, ey, ez ! reacion field |
82 |
double precision, dimension(natoms) :: mu ! dipole moments |
83 |
double precision, dimension(natoms) :: Axx, Axy, Axz ! rotational array |
84 |
double precision, dimension(natoms) :: Ayx, Ayy, Ayz ! |
85 |
double precision, dimension(natoms) :: Azx, Azy, Azz ! |
86 |
|
87 |
! local variables |
88 |
|
89 |
double precision rxi, ryi, rzi |
90 |
double precision rcutsq, rlstsq, rrfsq, rijsq, rij |
91 |
double precision rxij, ryij, rzij |
92 |
double precision prerf ! a reaction field preterm |
93 |
|
94 |
double precision, dimension(9,natoms) :: A |
95 |
|
96 |
integer :: i, j, k, ix, jx, nlist |
97 |
integer :: jbeg, jend, jnab |
98 |
|
99 |
logical :: exclude_temp1, exclude_temp2, exclude |
100 |
|
101 |
|
102 |
!******************************************************************* |
103 |
|
104 |
do i=1, natoms |
105 |
A(1,i) = Axx(i) |
106 |
A(2,i) = Axy(i) |
107 |
A(3,i) = Axz(i) |
108 |
|
109 |
A(4,i) = Ayx(i) |
110 |
A(5,i) = Ayy(i) |
111 |
A(6,i) = Ayz(i) |
112 |
|
113 |
A(7,i) = Azx(i) |
114 |
A(8,i) = Azy(i) |
115 |
A(9,i) = Azz(i) |
116 |
end do |
117 |
|
118 |
|
119 |
rlstsq = rlist * rlist |
120 |
rcutsq = rcut * rcut |
121 |
rrfsq = rrf * rrf |
122 |
|
123 |
prerf = 14.39257d0 * 2.0d0 * ( dielectric - 1.0d0 ) / & |
124 |
( ( 2.0d0 * dielectric + 1.0d0 ) * rrfsq * rrf ) |
125 |
|
126 |
! ** zero forces ** |
127 |
|
128 |
do I = 1, natoms |
129 |
|
130 |
fx(i) = 0.0d0 |
131 |
fy(i) = 0.0d0 |
132 |
fz(i) = 0.0d0 |
133 |
|
134 |
tx(i) = 0.0d0 |
135 |
ty(i) = 0.0d0 |
136 |
tz(i) = 0.0d0 |
137 |
|
138 |
ex(i) = 0.0d0 |
139 |
ey(i) = 0.0d0 |
140 |
ez(i) = 0.0d0 |
141 |
|
142 |
end do |
143 |
|
144 |
V = 0.0d0 |
145 |
|
146 |
if ( update ) then |
147 |
|
148 |
|
149 |
! ** save current configuration, construct ** |
150 |
! ** neighbour list and calculate forces ** |
151 |
|
152 |
call long_range_save(natoms,rx,ry,rz,rx0,ry0,rz0) |
153 |
|
154 |
nlist = 0 |
155 |
|
156 |
do i = 1, natoms - 1 |
157 |
|
158 |
point(i) = nlist + 1 |
159 |
|
160 |
rxi = rx(i) |
161 |
ryi = ry(i) |
162 |
rzi = rz(i) |
163 |
|
164 |
do j = i + 1, natoms |
165 |
|
166 |
rxij = rx(j) - rxi |
167 |
ryij = ry(j) - ryi |
168 |
rzij = rz(j) - rzi |
169 |
|
170 |
exclude = .false. |
171 |
|
172 |
if( check_exclude ) then |
173 |
|
174 |
do k = 1, n_pairs |
175 |
|
176 |
! add 1 for c -> fortran indexing |
177 |
ix = pair_i(k) + 1 |
178 |
jx = pair_j(k) + 1 |
179 |
exclude_temp1 = (i.eq.ix).and.(j.eq.jx) |
180 |
exclude_temp2 = (i.eq.jx).and.(j.eq.ix) |
181 |
|
182 |
if(exclude_temp1.or.exclude_temp2) then |
183 |
|
184 |
exclude = .true. |
185 |
endif |
186 |
enddo |
187 |
endif |
188 |
|
189 |
if(.not.exclude) then |
190 |
|
191 |
rxij = rxij - box_x * dsign( 1.0d0, rxij ) & |
192 |
* int( (dabs( rxij / box_x ) + 0.5d0) ) |
193 |
ryij = ryij - box_y * dsign( 1.0d0, ryij ) & |
194 |
* int( (dabs( ryij / box_y ) + 0.5d0) ) |
195 |
rzij = rzij - box_z * dsign( 1.0d0, rzij ) & |
196 |
* int( (dabs( rzij / box_z ) + 0.5d0) ) |
197 |
|
198 |
rijsq = rxij * rxij + ryij * ryij + rzij * rzij |
199 |
|
200 |
if ( rijsq .lt. rlstsq ) then |
201 |
|
202 |
nlist = nlist + 1 |
203 |
list(nlist) = j |
204 |
|
205 |
! ** remove this check if maxnab is appropriate ** |
206 |
|
207 |
if ( nlist .eq. maxnab ) then |
208 |
write(*,*) i, j, nlist, maxnab, rijsq, rlstsq |
209 |
stop 'list too small' |
210 |
endif |
211 |
|
212 |
if ( rijsq .lt. rcutsq ) then |
213 |
|
214 |
if ( is_vdw(i) .and. is_vdw(j) ) then |
215 |
call force_VDW( i, j, rcutsq, rijsq, sigma, epslon, & |
216 |
v, fx, fy, fz, rxij, ryij, rzij, natoms ) |
217 |
endif |
218 |
|
219 |
if ( is_lj(i) .and. is_lj(j) ) then |
220 |
call force_lj( i, j, rcutsq, rijsq, sigma, epslon, v, & |
221 |
fx, fy, fz, rxij, ryij, rzij, natoms ) |
222 |
endif |
223 |
|
224 |
if( is_ssd(i) .and. is_ssd(j) ) then |
225 |
|
226 |
rij = dsqrt(rijsq); |
227 |
|
228 |
call ssd_forces( natoms, i, j, i, j, rxij, ryij, rzij, & |
229 |
rij, v, rijsq, fx, fy, fz, tx, ty, tz, A ) |
230 |
end if |
231 |
|
232 |
endif |
233 |
|
234 |
if( is_dipole(i) .and. is_dipole(j) ) then |
235 |
if( rijsq .lt. rrfsq ) then |
236 |
|
237 |
call dipole_rf( i, j, rrfsq, rijsq, v, rxij, ryij, & |
238 |
rzij, rrf, dielectric, rtaper, natoms, & |
239 |
ux, uy, uz, mu, fx, fy, fz, tx, ty, tz, & |
240 |
ex, ey, ez, prerf ) |
241 |
endif |
242 |
end if |
243 |
endif |
244 |
endif |
245 |
|
246 |
enddo |
247 |
|
248 |
end do |
249 |
|
250 |
point(natoms) = nlist + 1 |
251 |
|
252 |
else |
253 |
|
254 |
|
255 |
!** use the list to find the neighbours ** |
256 |
|
257 |
do i = 1, natoms - 1 |
258 |
|
259 |
jbeg = point(i) |
260 |
jend = point(i+1) - 1 |
261 |
|
262 |
!** check that atom i has neighbours ** |
263 |
|
264 |
if( jbeg .le. jend ) then |
265 |
|
266 |
rxi = rx(i) |
267 |
ryi = ry(i) |
268 |
rzi = rz(i) |
269 |
|
270 |
do jnab = jbeg, jend |
271 |
|
272 |
j = list(jnab) |
273 |
|
274 |
rxij = rx(j) - rxi |
275 |
ryij = ry(j) - ryi |
276 |
rzij = rz(j) - rzi |
277 |
|
278 |
rxij = rxij - box_x * dsign( 1.0d0, rxij ) & |
279 |
* int( dabs( rxij / box_x ) + 0.5d0 ) |
280 |
ryij = ryij - box_y * dsign( 1.0d0, ryij ) & |
281 |
* int( dabs( ryij / box_y ) + 0.5d0 ) |
282 |
rzij = rzij - box_z * dsign( 1.0d0, rzij ) & |
283 |
* int( dabs( rzij / box_z ) + 0.5d0 ) |
284 |
|
285 |
rijsq = rxij * rxij + ryij * ryij + rzij * rzij |
286 |
|
287 |
if ( rijsq .lt. rcutsq ) then |
288 |
|
289 |
if ( is_vdw(i) .and. is_vdw(j) ) then |
290 |
call force_VDW( i, j, rcutsq, rijsq, sigma, epslon, & |
291 |
v, fx, fy, fz, rxij, ryij, rzij, natoms ) |
292 |
endif |
293 |
|
294 |
if ( is_lj(i) .and. is_lj(j) ) then |
295 |
call force_lj( i, j, rcutsq, rijsq, sigma, epslon, v, & |
296 |
fx, fy, fz, rxij, ryij, rzij, natoms ) |
297 |
endif |
298 |
|
299 |
if( is_ssd(i) .and. is_ssd(j) ) then |
300 |
|
301 |
rij = dsqrt(rijsq); |
302 |
|
303 |
call ssd_forces( natoms, i, j, i, j, rxij, ryij, rzij, & |
304 |
rij, v, rijsq, fx, fy, fz, tx, ty, tz, A ) |
305 |
end if |
306 |
|
307 |
endif |
308 |
|
309 |
|
310 |
if( is_dipole(i) .and. is_dipole(j) ) then |
311 |
if( rijsq .lt. rrfsq ) then |
312 |
|
313 |
call dipole_rf( i, j, rrfsq, rijsq, v, rxij, ryij, & |
314 |
rzij, rrf, dielectric, rtaper, natoms, & |
315 |
ux, uy, uz, mu, fx, fy, fz, tx, ty, tz, & |
316 |
ex, ey, ez, prerf ) |
317 |
endif |
318 |
endif |
319 |
|
320 |
end do |
321 |
endif |
322 |
end do |
323 |
endif |
324 |
|
325 |
! calculate the reaction field effect |
326 |
|
327 |
!!$ call reaction_field( v, prerf, natoms, mu, ux, uy, uz, tx, ty, tz, & |
328 |
!!$ ex, ey, ez, is_dipole ) |
329 |
|
330 |
|
331 |
end subroutine long_range_force |
332 |
|
333 |
|
334 |
|
335 |
subroutine long_range_check ( rcut, rlist, update, natoms, & |
336 |
rx,ry,rz,rx0,ry0,rz0) |
337 |
|
338 |
|
339 |
!******************************************************************* |
340 |
!** decides whether the list needs to be reconstructed. ** |
341 |
!** ** |
342 |
!** principal variables: ** |
343 |
!** ** |
344 |
!** real rx(n),ry(n),rz(n) atom positions ** |
345 |
!** real rx0(n),ry0(n),rz0(n) coordinates at last update ** |
346 |
!** real rlist radius of verlet list ** |
347 |
!** real rcut cutoff distance for forces ** |
348 |
!** real dispmx largest displacement ** |
349 |
!** integer n number of atoms ** |
350 |
!** logical update if true the list is updated ** |
351 |
!** ** |
352 |
!** usage: ** |
353 |
!** ** |
354 |
!** check is called to set update before every call to force. ** |
355 |
!******************************************************************* |
356 |
|
357 |
|
358 |
! Passed variables |
359 |
|
360 |
integer :: natoms |
361 |
|
362 |
logical(kind=2) :: update |
363 |
|
364 |
double precision rcut, rlist |
365 |
|
366 |
! Passed arrays |
367 |
double precision, dimension(natoms) :: rx, ry, rz |
368 |
double precision, dimension(natoms) :: rx0, ry0, rz0 |
369 |
|
370 |
! local variables |
371 |
|
372 |
double precision dispmx |
373 |
|
374 |
integer :: i |
375 |
|
376 |
!******************************************************************* |
377 |
|
378 |
!** calculate maximum displacement since last update ** |
379 |
|
380 |
dispmx = 0.0d0 |
381 |
|
382 |
do i = 1, natoms |
383 |
|
384 |
!write(*,*) 'calling displacement', i |
385 |
|
386 |
dispmx = max ( dabs ( rx(i) - rx0(i) ), dispmx ) |
387 |
dispmx = max ( dabs ( ry(i) - ry0(i) ), dispmx ) |
388 |
dispmx = max ( dabs ( rz(i) - rz0(i) ), dispmx ) |
389 |
|
390 |
!write(*,*) 'called displacement', i |
391 |
|
392 |
end do |
393 |
|
394 |
!** a conservative test of the list skin crossing ** |
395 |
|
396 |
dispmx = 2.0d0 * dsqrt ( 3.0d0 * dispmx ** 2 ) |
397 |
|
398 |
update = ( dispmx .gt. ( rlist - rcut ) ) |
399 |
|
400 |
|
401 |
return |
402 |
end subroutine long_range_check |
403 |
|
404 |
|
405 |
|
406 |
subroutine long_range_save(natoms,rx,ry,rz,rx0,ry0,rz0) |
407 |
|
408 |
|
409 |
!******************************************************************* |
410 |
!** saves current configuration for future checking. ** |
411 |
!** ** |
412 |
!** principal variables: ** |
413 |
!** ** |
414 |
!** real rx(n),ry(n),rz(n) atom positions ** |
415 |
!** real rx0(n),ry0(n),rz0(n) storage locations for save ** |
416 |
!** integer n number of atoms ** |
417 |
!** ** |
418 |
!** usage: ** |
419 |
!** ** |
420 |
!** save is called whenever the new verlet list is constructed. ** |
421 |
!******************************************************************* |
422 |
|
423 |
integer :: i |
424 |
|
425 |
integer :: natoms |
426 |
|
427 |
double precision, dimension(natoms) :: rx, ry, rz |
428 |
double precision, dimension(natoms) :: rx0, ry0, rz0 |
429 |
|
430 |
!******************************************************************* |
431 |
|
432 |
do i = 1, natoms |
433 |
|
434 |
rx0(i) = rx(i) |
435 |
ry0(i) = ry(i) |
436 |
rz0(i) = rz(i) |
437 |
|
438 |
end do |
439 |
|
440 |
return |
441 |
end subroutine long_range_save |
442 |
|
443 |
end module long_range_forces |