1 |
mmeineke |
377 |
!! do_Forces.F90 |
2 |
|
|
!! module do_Forces |
3 |
|
|
!! Calculates Long Range forces. |
4 |
|
|
|
5 |
|
|
!! @author Charles F. Vardeman II |
6 |
|
|
!! @author Matthew Meineke |
7 |
chuckv |
480 |
!! @version $Id: do_Forces.F90,v 1.12 2003-04-08 17:16:22 chuckv Exp $, $Date: 2003-04-08 17:16:22 $, $Name: not supported by cvs2svn $, $Revision: 1.12 $ |
8 |
mmeineke |
377 |
|
9 |
|
|
module do_Forces |
10 |
|
|
use force_globals |
11 |
|
|
use simulation |
12 |
|
|
use definitions |
13 |
|
|
use atype_module |
14 |
|
|
use neighborLists |
15 |
|
|
use lj |
16 |
|
|
use sticky_pair |
17 |
|
|
use dipole_dipole |
18 |
|
|
use reaction_field |
19 |
|
|
use gb_pair |
20 |
|
|
#ifdef IS_MPI |
21 |
|
|
use mpiSimulation |
22 |
|
|
#endif |
23 |
|
|
|
24 |
|
|
implicit none |
25 |
|
|
PRIVATE |
26 |
|
|
|
27 |
|
|
#define __FORTRAN90 |
28 |
|
|
#include "fForceField.h" |
29 |
|
|
|
30 |
|
|
logical, save :: do_forces_initialized = .false. |
31 |
|
|
logical, save :: FF_uses_LJ |
32 |
|
|
logical, save :: FF_uses_sticky |
33 |
|
|
logical, save :: FF_uses_dipoles |
34 |
|
|
logical, save :: FF_uses_RF |
35 |
|
|
logical, save :: FF_uses_GB |
36 |
|
|
logical, save :: FF_uses_EAM |
37 |
|
|
|
38 |
|
|
public :: init_FF |
39 |
|
|
public :: do_force_loop |
40 |
|
|
|
41 |
|
|
contains |
42 |
|
|
|
43 |
|
|
subroutine init_FF(LJMIXPOLICY, use_RF_c, thisStat) |
44 |
|
|
|
45 |
|
|
integer, intent(in) :: LJMIXPOLICY |
46 |
|
|
logical, intent(in) :: use_RF_c |
47 |
|
|
|
48 |
|
|
integer, intent(out) :: thisStat |
49 |
|
|
integer :: my_status, nMatches |
50 |
|
|
integer, pointer :: MatchList(:) => null() |
51 |
|
|
real(kind=dp) :: rcut, rrf, rt, dielect |
52 |
|
|
|
53 |
|
|
!! assume things are copacetic, unless they aren't |
54 |
|
|
thisStat = 0 |
55 |
|
|
|
56 |
|
|
!! Fortran's version of a cast: |
57 |
|
|
FF_uses_RF = use_RF_c |
58 |
|
|
|
59 |
|
|
!! init_FF is called *after* all of the atom types have been |
60 |
|
|
!! defined in atype_module using the new_atype subroutine. |
61 |
|
|
!! |
62 |
|
|
!! this will scan through the known atypes and figure out what |
63 |
|
|
!! interactions are used by the force field. |
64 |
|
|
|
65 |
|
|
FF_uses_LJ = .false. |
66 |
|
|
FF_uses_sticky = .false. |
67 |
|
|
FF_uses_dipoles = .false. |
68 |
|
|
FF_uses_GB = .false. |
69 |
|
|
FF_uses_EAM = .false. |
70 |
|
|
|
71 |
|
|
call getMatchingElementList(atypes, "is_LJ", .true., nMatches, MatchList) |
72 |
|
|
if (nMatches .gt. 0) FF_uses_LJ = .true. |
73 |
|
|
|
74 |
|
|
call getMatchingElementList(atypes, "is_DP", .true., nMatches, MatchList) |
75 |
|
|
if (nMatches .gt. 0) FF_uses_dipoles = .true. |
76 |
|
|
|
77 |
|
|
call getMatchingElementList(atypes, "is_Sticky", .true., nMatches, & |
78 |
|
|
MatchList) |
79 |
|
|
if (nMatches .gt. 0) FF_uses_Sticky = .true. |
80 |
|
|
|
81 |
|
|
call getMatchingElementList(atypes, "is_GB", .true., nMatches, MatchList) |
82 |
|
|
if (nMatches .gt. 0) FF_uses_GB = .true. |
83 |
|
|
|
84 |
|
|
call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList) |
85 |
|
|
if (nMatches .gt. 0) FF_uses_EAM = .true. |
86 |
|
|
|
87 |
|
|
!! check to make sure the FF_uses_RF setting makes sense |
88 |
|
|
|
89 |
|
|
if (FF_uses_dipoles) then |
90 |
|
|
rrf = getRrf() |
91 |
|
|
rt = getRt() |
92 |
|
|
call initialize_dipole(rrf, rt) |
93 |
|
|
if (FF_uses_RF) then |
94 |
|
|
dielect = getDielect() |
95 |
|
|
call initialize_rf(rrf, rt, dielect) |
96 |
|
|
endif |
97 |
|
|
else |
98 |
|
|
if (FF_uses_RF) then |
99 |
|
|
write(default_error,*) 'Using Reaction Field with no dipoles? Huh?' |
100 |
|
|
thisStat = -1 |
101 |
|
|
return |
102 |
|
|
endif |
103 |
|
|
endif |
104 |
|
|
|
105 |
|
|
if (FF_uses_LJ) then |
106 |
|
|
|
107 |
|
|
call getRcut(rcut) |
108 |
|
|
|
109 |
|
|
select case (LJMIXPOLICY) |
110 |
|
|
case (LB_MIXING_RULE) |
111 |
|
|
call init_lj_FF(LB_MIXING_RULE, rcut, my_status) |
112 |
|
|
case (EXPLICIT_MIXING_RULE) |
113 |
|
|
call init_lj_FF(EXPLICIT_MIXING_RULE, rcut, my_status) |
114 |
|
|
case default |
115 |
|
|
write(default_error,*) 'unknown LJ Mixing Policy!' |
116 |
|
|
thisStat = -1 |
117 |
|
|
return |
118 |
|
|
end select |
119 |
|
|
if (my_status /= 0) then |
120 |
|
|
thisStat = -1 |
121 |
|
|
return |
122 |
|
|
end if |
123 |
|
|
endif |
124 |
|
|
|
125 |
|
|
if (FF_uses_sticky) then |
126 |
|
|
call check_sticky_FF(my_status) |
127 |
|
|
if (my_status /= 0) then |
128 |
|
|
thisStat = -1 |
129 |
|
|
return |
130 |
|
|
end if |
131 |
|
|
endif |
132 |
|
|
|
133 |
|
|
if (FF_uses_GB) then |
134 |
|
|
call check_gb_pair_FF(my_status) |
135 |
|
|
if (my_status .ne. 0) then |
136 |
|
|
thisStat = -1 |
137 |
|
|
return |
138 |
|
|
endif |
139 |
|
|
endif |
140 |
|
|
|
141 |
|
|
if (FF_uses_GB .and. FF_uses_LJ) then |
142 |
|
|
endif |
143 |
chuckv |
480 |
if (.not. do_forces_initialized) then |
144 |
|
|
!! Create neighbor lists |
145 |
|
|
call expandNeighborList(getNlocal(), my_status) |
146 |
|
|
if (my_Status /= 0) then |
147 |
|
|
write(default_error,*) "SimSetup: ExpandNeighborList returned error." |
148 |
|
|
thisStat = -1 |
149 |
|
|
return |
150 |
|
|
endif |
151 |
|
|
endif |
152 |
mmeineke |
377 |
|
153 |
|
|
do_forces_initialized = .true. |
154 |
|
|
|
155 |
|
|
end subroutine init_FF |
156 |
|
|
|
157 |
|
|
|
158 |
|
|
!! Does force loop over i,j pairs. Calls do_pair to calculates forces. |
159 |
|
|
!-------------------------------------------------------------> |
160 |
|
|
subroutine do_force_loop(q, A, u_l, f, t, tau, pot, do_pot_c, do_stress_c, & |
161 |
|
|
error) |
162 |
|
|
!! Position array provided by C, dimensioned by getNlocal |
163 |
|
|
real ( kind = dp ), dimension(3,getNlocal()) :: q |
164 |
|
|
!! Rotation Matrix for each long range particle in simulation. |
165 |
|
|
real( kind = dp), dimension(9,getNlocal()) :: A |
166 |
|
|
!! Unit vectors for dipoles (lab frame) |
167 |
|
|
real( kind = dp ), dimension(3,getNlocal()) :: u_l |
168 |
|
|
!! Force array provided by C, dimensioned by getNlocal |
169 |
|
|
real ( kind = dp ), dimension(3,getNlocal()) :: f |
170 |
|
|
!! Torsion array provided by C, dimensioned by getNlocal |
171 |
|
|
real( kind = dp ), dimension(3,getNlocal()) :: t |
172 |
|
|
!! Stress Tensor |
173 |
|
|
real( kind = dp), dimension(9) :: tau |
174 |
|
|
real ( kind = dp ) :: pot |
175 |
|
|
logical ( kind = 2) :: do_pot_c, do_stress_c |
176 |
|
|
logical :: do_pot |
177 |
|
|
logical :: do_stress |
178 |
chuckv |
439 |
#ifdef IS_MPI |
179 |
chuckv |
441 |
real( kind = DP ) :: pot_local |
180 |
mmeineke |
377 |
integer :: nrow |
181 |
|
|
integer :: ncol |
182 |
|
|
#endif |
183 |
|
|
integer :: nlocal |
184 |
|
|
integer :: natoms |
185 |
|
|
logical :: update_nlist |
186 |
|
|
integer :: i, j, jbeg, jend, jnab |
187 |
|
|
integer :: nlist |
188 |
|
|
real( kind = DP ) :: rijsq, rlistsq, rcutsq, rlist, rcut |
189 |
|
|
real(kind=dp),dimension(3) :: d |
190 |
|
|
real(kind=dp) :: rfpot, mu_i, virial |
191 |
|
|
integer :: me_i |
192 |
|
|
logical :: is_dp_i |
193 |
|
|
integer :: neighborListSize |
194 |
|
|
integer :: listerror, error |
195 |
|
|
integer :: localError |
196 |
|
|
|
197 |
|
|
!! initialize local variables |
198 |
|
|
|
199 |
|
|
#ifdef IS_MPI |
200 |
chuckv |
441 |
pot_local = 0.0_dp |
201 |
mmeineke |
377 |
nlocal = getNlocal() |
202 |
|
|
nrow = getNrow(plan_row) |
203 |
|
|
ncol = getNcol(plan_col) |
204 |
|
|
#else |
205 |
|
|
nlocal = getNlocal() |
206 |
|
|
natoms = nlocal |
207 |
|
|
#endif |
208 |
chuckv |
441 |
|
209 |
mmeineke |
377 |
call getRcut(rcut,rc2=rcutsq) |
210 |
|
|
call getRlist(rlist,rlistsq) |
211 |
|
|
|
212 |
|
|
call check_initialization(localError) |
213 |
|
|
if ( localError .ne. 0 ) then |
214 |
|
|
error = -1 |
215 |
|
|
return |
216 |
|
|
end if |
217 |
|
|
call zero_work_arrays() |
218 |
|
|
|
219 |
|
|
do_pot = do_pot_c |
220 |
|
|
do_stress = do_stress_c |
221 |
chuckv |
470 |
|
222 |
mmeineke |
377 |
|
223 |
|
|
! Gather all information needed by all force loops: |
224 |
|
|
|
225 |
|
|
#ifdef IS_MPI |
226 |
|
|
|
227 |
|
|
call gather(q,q_Row,plan_row3d) |
228 |
|
|
call gather(q,q_Col,plan_col3d) |
229 |
|
|
|
230 |
|
|
if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then |
231 |
|
|
call gather(u_l,u_l_Row,plan_row3d) |
232 |
|
|
call gather(u_l,u_l_Col,plan_col3d) |
233 |
|
|
|
234 |
|
|
call gather(A,A_Row,plan_row_rotation) |
235 |
|
|
call gather(A,A_Col,plan_col_rotation) |
236 |
|
|
endif |
237 |
|
|
|
238 |
|
|
#endif |
239 |
|
|
|
240 |
|
|
if (FF_RequiresPrepairCalc() .and. SimRequiresPrepairCalc()) then |
241 |
|
|
!! See if we need to update neighbor lists |
242 |
|
|
call checkNeighborList(nlocal, q, rcut, rlist, update_nlist) |
243 |
|
|
!! if_mpi_gather_stuff_for_prepair |
244 |
|
|
!! do_prepair_loop_if_needed |
245 |
|
|
!! if_mpi_scatter_stuff_from_prepair |
246 |
|
|
!! if_mpi_gather_stuff_from_prepair_to_main_loop |
247 |
|
|
else |
248 |
|
|
!! See if we need to update neighbor lists |
249 |
|
|
call checkNeighborList(nlocal, q, rcut, rlist, update_nlist) |
250 |
|
|
endif |
251 |
|
|
|
252 |
|
|
#ifdef IS_MPI |
253 |
|
|
|
254 |
|
|
if (update_nlist) then |
255 |
|
|
|
256 |
|
|
!! save current configuration, construct neighbor list, |
257 |
|
|
!! and calculate forces |
258 |
mmeineke |
459 |
call saveNeighborList(nlocal, q) |
259 |
mmeineke |
377 |
|
260 |
|
|
neighborListSize = size(list) |
261 |
|
|
nlist = 0 |
262 |
|
|
|
263 |
|
|
do i = 1, nrow |
264 |
|
|
point(i) = nlist + 1 |
265 |
|
|
|
266 |
|
|
inner: do j = 1, ncol |
267 |
|
|
|
268 |
|
|
if (skipThisPair(i,j)) cycle inner |
269 |
|
|
|
270 |
|
|
call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq) |
271 |
|
|
|
272 |
|
|
if (rijsq < rlistsq) then |
273 |
|
|
|
274 |
|
|
nlist = nlist + 1 |
275 |
|
|
|
276 |
|
|
if (nlist > neighborListSize) then |
277 |
|
|
call expandNeighborList(nlocal, listerror) |
278 |
|
|
if (listerror /= 0) then |
279 |
|
|
error = -1 |
280 |
|
|
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
281 |
|
|
return |
282 |
|
|
end if |
283 |
|
|
neighborListSize = size(list) |
284 |
|
|
endif |
285 |
|
|
|
286 |
|
|
list(nlist) = j |
287 |
|
|
|
288 |
|
|
if (rijsq < rcutsq) then |
289 |
|
|
call do_pair(i, j, rijsq, d, do_pot, do_stress, & |
290 |
chuckv |
441 |
u_l, A, f, t, pot_local) |
291 |
mmeineke |
377 |
endif |
292 |
|
|
endif |
293 |
|
|
enddo inner |
294 |
|
|
enddo |
295 |
|
|
|
296 |
|
|
point(nrow + 1) = nlist + 1 |
297 |
|
|
|
298 |
|
|
else !! (of update_check) |
299 |
|
|
|
300 |
|
|
! use the list to find the neighbors |
301 |
|
|
do i = 1, nrow |
302 |
|
|
JBEG = POINT(i) |
303 |
|
|
JEND = POINT(i+1) - 1 |
304 |
|
|
! check thiat molecule i has neighbors |
305 |
|
|
if (jbeg .le. jend) then |
306 |
|
|
|
307 |
|
|
do jnab = jbeg, jend |
308 |
|
|
j = list(jnab) |
309 |
|
|
|
310 |
|
|
call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq) |
311 |
|
|
call do_pair(i, j, rijsq, d, do_pot, do_stress, & |
312 |
chuckv |
441 |
u_l, A, f, t, pot_local) |
313 |
mmeineke |
377 |
|
314 |
|
|
enddo |
315 |
|
|
endif |
316 |
|
|
enddo |
317 |
|
|
endif |
318 |
|
|
|
319 |
|
|
#else |
320 |
|
|
|
321 |
|
|
if (update_nlist) then |
322 |
|
|
|
323 |
|
|
! save current configuration, contruct neighbor list, |
324 |
|
|
! and calculate forces |
325 |
mmeineke |
459 |
call saveNeighborList(natoms, q) |
326 |
mmeineke |
377 |
|
327 |
|
|
neighborListSize = size(list) |
328 |
|
|
|
329 |
|
|
nlist = 0 |
330 |
|
|
|
331 |
|
|
do i = 1, natoms-1 |
332 |
|
|
point(i) = nlist + 1 |
333 |
|
|
|
334 |
|
|
inner: do j = i+1, natoms |
335 |
|
|
|
336 |
chuckv |
388 |
if (skipThisPair(i,j)) cycle inner |
337 |
|
|
|
338 |
mmeineke |
377 |
call get_interatomic_vector(q(:,i), q(:,j), d, rijsq) |
339 |
|
|
|
340 |
chuckv |
388 |
|
341 |
mmeineke |
377 |
if (rijsq < rlistsq) then |
342 |
|
|
|
343 |
|
|
nlist = nlist + 1 |
344 |
|
|
|
345 |
|
|
if (nlist > neighborListSize) then |
346 |
|
|
call expandNeighborList(natoms, listerror) |
347 |
|
|
if (listerror /= 0) then |
348 |
|
|
error = -1 |
349 |
|
|
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
350 |
|
|
return |
351 |
|
|
end if |
352 |
|
|
neighborListSize = size(list) |
353 |
|
|
endif |
354 |
|
|
|
355 |
|
|
list(nlist) = j |
356 |
|
|
|
357 |
|
|
if (rijsq < rcutsq) then |
358 |
|
|
call do_pair(i, j, rijsq, d, do_pot, do_stress, & |
359 |
chuckv |
441 |
u_l, A, f, t, pot) |
360 |
mmeineke |
377 |
endif |
361 |
|
|
endif |
362 |
|
|
enddo inner |
363 |
|
|
enddo |
364 |
|
|
|
365 |
|
|
point(natoms) = nlist + 1 |
366 |
|
|
|
367 |
|
|
else !! (update) |
368 |
|
|
|
369 |
|
|
! use the list to find the neighbors |
370 |
|
|
do i = 1, natoms-1 |
371 |
|
|
JBEG = POINT(i) |
372 |
|
|
JEND = POINT(i+1) - 1 |
373 |
|
|
! check thiat molecule i has neighbors |
374 |
|
|
if (jbeg .le. jend) then |
375 |
|
|
|
376 |
|
|
do jnab = jbeg, jend |
377 |
|
|
j = list(jnab) |
378 |
|
|
|
379 |
|
|
call get_interatomic_vector(q(:,i), q(:,j), d, rijsq) |
380 |
|
|
call do_pair(i, j, rijsq, d, do_pot, do_stress, & |
381 |
chuckv |
441 |
u_l, A, f, t, pot) |
382 |
mmeineke |
377 |
|
383 |
|
|
enddo |
384 |
|
|
endif |
385 |
|
|
enddo |
386 |
|
|
endif |
387 |
|
|
|
388 |
|
|
#endif |
389 |
|
|
|
390 |
|
|
! phew, done with main loop. |
391 |
|
|
|
392 |
|
|
#ifdef IS_MPI |
393 |
|
|
!!distribute forces |
394 |
chuckv |
438 |
|
395 |
|
|
f_temp = 0.0_dp |
396 |
|
|
call scatter(f_Row,f_temp,plan_row3d) |
397 |
|
|
do i = 1,nlocal |
398 |
|
|
f(1:3,i) = f(1:3,i) + f_temp(1:3,i) |
399 |
|
|
end do |
400 |
|
|
|
401 |
|
|
f_temp = 0.0_dp |
402 |
mmeineke |
377 |
call scatter(f_Col,f_temp,plan_col3d) |
403 |
|
|
do i = 1,nlocal |
404 |
|
|
f(1:3,i) = f(1:3,i) + f_temp(1:3,i) |
405 |
|
|
end do |
406 |
|
|
|
407 |
|
|
if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then |
408 |
chuckv |
438 |
t_temp = 0.0_dp |
409 |
|
|
call scatter(t_Row,t_temp,plan_row3d) |
410 |
|
|
do i = 1,nlocal |
411 |
|
|
t(1:3,i) = t(1:3,i) + t_temp(1:3,i) |
412 |
|
|
end do |
413 |
|
|
t_temp = 0.0_dp |
414 |
mmeineke |
377 |
call scatter(t_Col,t_temp,plan_col3d) |
415 |
|
|
|
416 |
|
|
do i = 1,nlocal |
417 |
|
|
t(1:3,i) = t(1:3,i) + t_temp(1:3,i) |
418 |
|
|
end do |
419 |
|
|
endif |
420 |
|
|
|
421 |
|
|
if (do_pot) then |
422 |
|
|
! scatter/gather pot_row into the members of my column |
423 |
|
|
call scatter(pot_Row, pot_Temp, plan_row) |
424 |
chuckv |
439 |
|
425 |
mmeineke |
377 |
! scatter/gather pot_local into all other procs |
426 |
|
|
! add resultant to get total pot |
427 |
|
|
do i = 1, nlocal |
428 |
|
|
pot_local = pot_local + pot_Temp(i) |
429 |
|
|
enddo |
430 |
chuckv |
439 |
|
431 |
|
|
pot_Temp = 0.0_DP |
432 |
mmeineke |
377 |
|
433 |
|
|
call scatter(pot_Col, pot_Temp, plan_col) |
434 |
|
|
do i = 1, nlocal |
435 |
|
|
pot_local = pot_local + pot_Temp(i) |
436 |
|
|
enddo |
437 |
chuckv |
439 |
|
438 |
mmeineke |
377 |
endif |
439 |
|
|
#endif |
440 |
|
|
|
441 |
|
|
if (FF_RequiresPostpairCalc() .and. SimRequiresPostpairCalc()) then |
442 |
|
|
|
443 |
|
|
if (FF_uses_RF .and. SimUsesRF()) then |
444 |
|
|
|
445 |
|
|
#ifdef IS_MPI |
446 |
|
|
call scatter(rf_Row,rf,plan_row3d) |
447 |
|
|
call scatter(rf_Col,rf_Temp,plan_col3d) |
448 |
|
|
do i = 1,nlocal |
449 |
|
|
rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i) |
450 |
|
|
end do |
451 |
|
|
#endif |
452 |
|
|
|
453 |
|
|
do i = 1, getNlocal() |
454 |
|
|
|
455 |
|
|
rfpot = 0.0_DP |
456 |
|
|
#ifdef IS_MPI |
457 |
|
|
me_i = atid_row(i) |
458 |
|
|
#else |
459 |
|
|
me_i = atid(i) |
460 |
|
|
#endif |
461 |
|
|
call getElementProperty(atypes, me_i, "is_DP", is_DP_i) |
462 |
|
|
if ( is_DP_i ) then |
463 |
|
|
call getElementProperty(atypes, me_i, "dipole_moment", mu_i) |
464 |
|
|
!! The reaction field needs to include a self contribution |
465 |
|
|
!! to the field: |
466 |
|
|
call accumulate_self_rf(i, mu_i, u_l) |
467 |
|
|
!! Get the reaction field contribution to the |
468 |
|
|
!! potential and torques: |
469 |
|
|
call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot) |
470 |
|
|
#ifdef IS_MPI |
471 |
|
|
pot_local = pot_local + rfpot |
472 |
|
|
#else |
473 |
|
|
pot = pot + rfpot |
474 |
|
|
|
475 |
|
|
#endif |
476 |
|
|
endif |
477 |
|
|
enddo |
478 |
|
|
endif |
479 |
|
|
endif |
480 |
|
|
|
481 |
|
|
|
482 |
|
|
#ifdef IS_MPI |
483 |
|
|
|
484 |
|
|
if (do_pot) then |
485 |
chuckv |
441 |
pot = pot + pot_local |
486 |
mmeineke |
377 |
!! we assume the c code will do the allreduce to get the total potential |
487 |
|
|
!! we could do it right here if we needed to... |
488 |
|
|
endif |
489 |
|
|
|
490 |
|
|
if (do_stress) then |
491 |
chuckv |
470 |
call mpi_allreduce(tau_Temp, tau,9,mpi_double_precision,mpi_sum, & |
492 |
mmeineke |
377 |
mpi_comm_world,mpi_err) |
493 |
chuckv |
470 |
call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, & |
494 |
mmeineke |
377 |
mpi_comm_world,mpi_err) |
495 |
|
|
endif |
496 |
|
|
|
497 |
|
|
#else |
498 |
|
|
|
499 |
|
|
if (do_stress) then |
500 |
|
|
tau = tau_Temp |
501 |
|
|
virial = virial_Temp |
502 |
|
|
endif |
503 |
|
|
|
504 |
|
|
#endif |
505 |
|
|
|
506 |
|
|
end subroutine do_force_loop |
507 |
|
|
|
508 |
chuckv |
441 |
subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot) |
509 |
mmeineke |
377 |
|
510 |
|
|
real( kind = dp ) :: pot |
511 |
chuckv |
460 |
real( kind = dp ), dimension(3,getNlocal()) :: u_l |
512 |
|
|
real (kind=dp), dimension(9,getNlocal()) :: A |
513 |
|
|
real (kind=dp), dimension(3,getNlocal()) :: f |
514 |
|
|
real (kind=dp), dimension(3,getNlocal()) :: t |
515 |
mmeineke |
377 |
|
516 |
|
|
logical, intent(inout) :: do_pot, do_stress |
517 |
|
|
integer, intent(in) :: i, j |
518 |
|
|
real ( kind = dp ), intent(inout) :: rijsq |
519 |
|
|
real ( kind = dp ) :: r |
520 |
|
|
real ( kind = dp ), intent(inout) :: d(3) |
521 |
|
|
logical :: is_LJ_i, is_LJ_j |
522 |
|
|
logical :: is_DP_i, is_DP_j |
523 |
|
|
logical :: is_GB_i, is_GB_j |
524 |
|
|
logical :: is_Sticky_i, is_Sticky_j |
525 |
|
|
integer :: me_i, me_j |
526 |
|
|
|
527 |
|
|
r = sqrt(rijsq) |
528 |
|
|
|
529 |
|
|
#ifdef IS_MPI |
530 |
|
|
|
531 |
|
|
me_i = atid_row(i) |
532 |
|
|
me_j = atid_col(j) |
533 |
|
|
|
534 |
|
|
#else |
535 |
|
|
|
536 |
|
|
me_i = atid(i) |
537 |
|
|
me_j = atid(j) |
538 |
|
|
|
539 |
|
|
#endif |
540 |
|
|
|
541 |
|
|
if (FF_uses_LJ .and. SimUsesLJ()) then |
542 |
|
|
call getElementProperty(atypes, me_i, "is_LJ", is_LJ_i) |
543 |
|
|
call getElementProperty(atypes, me_j, "is_LJ", is_LJ_j) |
544 |
|
|
|
545 |
|
|
if ( is_LJ_i .and. is_LJ_j ) & |
546 |
|
|
call do_lj_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress) |
547 |
|
|
endif |
548 |
|
|
|
549 |
|
|
if (FF_uses_dipoles .and. SimUsesDipoles()) then |
550 |
|
|
call getElementProperty(atypes, me_i, "is_DP", is_DP_i) |
551 |
|
|
call getElementProperty(atypes, me_j, "is_DP", is_DP_j) |
552 |
|
|
|
553 |
|
|
if ( is_DP_i .and. is_DP_j ) then |
554 |
|
|
|
555 |
gezelter |
462 |
call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, & |
556 |
mmeineke |
377 |
do_pot, do_stress) |
557 |
|
|
if (FF_uses_RF .and. SimUsesRF()) then |
558 |
|
|
call accumulate_rf(i, j, r, u_l) |
559 |
|
|
call rf_correct_forces(i, j, d, r, u_l, f, do_stress) |
560 |
|
|
endif |
561 |
|
|
|
562 |
|
|
endif |
563 |
|
|
endif |
564 |
|
|
|
565 |
|
|
if (FF_uses_Sticky .and. SimUsesSticky()) then |
566 |
|
|
|
567 |
|
|
call getElementProperty(atypes, me_i, "is_Sticky", is_Sticky_i) |
568 |
|
|
call getElementProperty(atypes, me_j, "is_Sticky", is_Sticky_j) |
569 |
chuckv |
388 |
|
570 |
mmeineke |
377 |
if ( is_Sticky_i .and. is_Sticky_j ) then |
571 |
|
|
call do_sticky_pair(i, j, d, r, rijsq, A, pot, f, t, & |
572 |
|
|
do_pot, do_stress) |
573 |
|
|
endif |
574 |
|
|
endif |
575 |
|
|
|
576 |
|
|
|
577 |
|
|
if (FF_uses_GB .and. SimUsesGB()) then |
578 |
|
|
|
579 |
|
|
call getElementProperty(atypes, me_i, "is_GB", is_GB_i) |
580 |
|
|
call getElementProperty(atypes, me_j, "is_GB", is_GB_j) |
581 |
|
|
|
582 |
|
|
if ( is_GB_i .and. is_GB_j ) then |
583 |
|
|
call do_gb_pair(i, j, d, r, rijsq, u_l, pot, f, t, & |
584 |
|
|
do_pot, do_stress) |
585 |
|
|
endif |
586 |
|
|
endif |
587 |
|
|
|
588 |
|
|
end subroutine do_pair |
589 |
|
|
|
590 |
|
|
|
591 |
|
|
subroutine get_interatomic_vector(q_i, q_j, d, r_sq) |
592 |
|
|
|
593 |
|
|
real (kind = dp), dimension(3) :: q_i |
594 |
|
|
real (kind = dp), dimension(3) :: q_j |
595 |
|
|
real ( kind = dp ), intent(out) :: r_sq |
596 |
|
|
real( kind = dp ) :: d(3) |
597 |
|
|
real( kind = dp ) :: d_old(3) |
598 |
|
|
d(1:3) = q_i(1:3) - q_j(1:3) |
599 |
|
|
d_old = d |
600 |
|
|
! Wrap back into periodic box if necessary |
601 |
|
|
if ( SimUsesPBC() ) then |
602 |
mmeineke |
393 |
|
603 |
mmeineke |
377 |
d(1:3) = d(1:3) - box(1:3) * sign(1.0_dp,d(1:3)) * & |
604 |
|
|
int(abs(d(1:3)/box(1:3)) + 0.5_dp) |
605 |
mmeineke |
393 |
|
606 |
mmeineke |
377 |
endif |
607 |
|
|
r_sq = dot_product(d,d) |
608 |
|
|
|
609 |
|
|
end subroutine get_interatomic_vector |
610 |
|
|
|
611 |
|
|
subroutine check_initialization(error) |
612 |
|
|
integer, intent(out) :: error |
613 |
|
|
|
614 |
|
|
error = 0 |
615 |
|
|
! Make sure we are properly initialized. |
616 |
|
|
if (.not. do_forces_initialized) then |
617 |
|
|
error = -1 |
618 |
|
|
return |
619 |
|
|
endif |
620 |
|
|
|
621 |
|
|
#ifdef IS_MPI |
622 |
|
|
if (.not. isMPISimSet()) then |
623 |
|
|
write(default_error,*) "ERROR: mpiSimulation has not been initialized!" |
624 |
|
|
error = -1 |
625 |
|
|
return |
626 |
|
|
endif |
627 |
|
|
#endif |
628 |
|
|
|
629 |
|
|
return |
630 |
|
|
end subroutine check_initialization |
631 |
|
|
|
632 |
|
|
|
633 |
|
|
subroutine zero_work_arrays() |
634 |
|
|
|
635 |
|
|
#ifdef IS_MPI |
636 |
|
|
|
637 |
|
|
q_Row = 0.0_dp |
638 |
|
|
q_Col = 0.0_dp |
639 |
|
|
|
640 |
|
|
u_l_Row = 0.0_dp |
641 |
|
|
u_l_Col = 0.0_dp |
642 |
|
|
|
643 |
|
|
A_Row = 0.0_dp |
644 |
|
|
A_Col = 0.0_dp |
645 |
|
|
|
646 |
|
|
f_Row = 0.0_dp |
647 |
|
|
f_Col = 0.0_dp |
648 |
|
|
f_Temp = 0.0_dp |
649 |
|
|
|
650 |
|
|
t_Row = 0.0_dp |
651 |
|
|
t_Col = 0.0_dp |
652 |
|
|
t_Temp = 0.0_dp |
653 |
|
|
|
654 |
|
|
pot_Row = 0.0_dp |
655 |
|
|
pot_Col = 0.0_dp |
656 |
|
|
pot_Temp = 0.0_dp |
657 |
|
|
|
658 |
|
|
rf_Row = 0.0_dp |
659 |
|
|
rf_Col = 0.0_dp |
660 |
|
|
rf_Temp = 0.0_dp |
661 |
|
|
|
662 |
|
|
#endif |
663 |
|
|
|
664 |
|
|
rf = 0.0_dp |
665 |
|
|
tau_Temp = 0.0_dp |
666 |
|
|
virial_Temp = 0.0_dp |
667 |
|
|
end subroutine zero_work_arrays |
668 |
|
|
|
669 |
|
|
function skipThisPair(atom1, atom2) result(skip_it) |
670 |
|
|
integer, intent(in) :: atom1 |
671 |
|
|
integer, intent(in), optional :: atom2 |
672 |
|
|
logical :: skip_it |
673 |
|
|
integer :: unique_id_1, unique_id_2 |
674 |
chuckv |
388 |
integer :: me_i,me_j |
675 |
mmeineke |
377 |
integer :: i |
676 |
|
|
|
677 |
|
|
skip_it = .false. |
678 |
|
|
|
679 |
|
|
!! there are a number of reasons to skip a pair or a particle |
680 |
|
|
!! mostly we do this to exclude atoms who are involved in short |
681 |
|
|
!! range interactions (bonds, bends, torsions), but we also need |
682 |
|
|
!! to exclude some overcounted interactions that result from |
683 |
|
|
!! the parallel decomposition |
684 |
|
|
|
685 |
|
|
#ifdef IS_MPI |
686 |
|
|
!! in MPI, we have to look up the unique IDs for each atom |
687 |
|
|
unique_id_1 = tagRow(atom1) |
688 |
|
|
#else |
689 |
|
|
!! in the normal loop, the atom numbers are unique |
690 |
|
|
unique_id_1 = atom1 |
691 |
|
|
#endif |
692 |
chuckv |
388 |
|
693 |
mmeineke |
377 |
!! We were called with only one atom, so just check the global exclude |
694 |
|
|
!! list for this atom |
695 |
|
|
if (.not. present(atom2)) then |
696 |
|
|
do i = 1, nExcludes_global |
697 |
|
|
if (excludesGlobal(i) == unique_id_1) then |
698 |
|
|
skip_it = .true. |
699 |
|
|
return |
700 |
|
|
end if |
701 |
|
|
end do |
702 |
|
|
return |
703 |
|
|
end if |
704 |
|
|
|
705 |
|
|
#ifdef IS_MPI |
706 |
|
|
unique_id_2 = tagColumn(atom2) |
707 |
|
|
#else |
708 |
|
|
unique_id_2 = atom2 |
709 |
|
|
#endif |
710 |
chuckv |
441 |
|
711 |
mmeineke |
377 |
#ifdef IS_MPI |
712 |
|
|
!! this situation should only arise in MPI simulations |
713 |
|
|
if (unique_id_1 == unique_id_2) then |
714 |
|
|
skip_it = .true. |
715 |
|
|
return |
716 |
|
|
end if |
717 |
|
|
|
718 |
|
|
!! this prevents us from doing the pair on multiple processors |
719 |
|
|
if (unique_id_1 < unique_id_2) then |
720 |
chuckv |
441 |
if (mod(unique_id_1 + unique_id_2,2) == 0) then |
721 |
|
|
skip_it = .true. |
722 |
|
|
return |
723 |
|
|
endif |
724 |
mmeineke |
377 |
else |
725 |
chuckv |
441 |
if (mod(unique_id_1 + unique_id_2,2) == 1) then |
726 |
|
|
skip_it = .true. |
727 |
|
|
return |
728 |
|
|
endif |
729 |
mmeineke |
377 |
endif |
730 |
|
|
#endif |
731 |
chuckv |
441 |
|
732 |
mmeineke |
377 |
!! the rest of these situations can happen in all simulations: |
733 |
|
|
do i = 1, nExcludes_global |
734 |
|
|
if ((excludesGlobal(i) == unique_id_1) .or. & |
735 |
|
|
(excludesGlobal(i) == unique_id_2)) then |
736 |
|
|
skip_it = .true. |
737 |
|
|
return |
738 |
|
|
endif |
739 |
|
|
enddo |
740 |
chuckv |
441 |
|
741 |
mmeineke |
377 |
do i = 1, nExcludes_local |
742 |
|
|
if (excludesLocal(1,i) == unique_id_1) then |
743 |
|
|
if (excludesLocal(2,i) == unique_id_2) then |
744 |
|
|
skip_it = .true. |
745 |
|
|
return |
746 |
|
|
endif |
747 |
|
|
else |
748 |
|
|
if (excludesLocal(1,i) == unique_id_2) then |
749 |
|
|
if (excludesLocal(2,i) == unique_id_1) then |
750 |
|
|
skip_it = .true. |
751 |
|
|
return |
752 |
|
|
endif |
753 |
|
|
endif |
754 |
|
|
endif |
755 |
|
|
end do |
756 |
|
|
|
757 |
|
|
return |
758 |
|
|
end function skipThisPair |
759 |
|
|
|
760 |
|
|
function FF_UsesDirectionalAtoms() result(doesit) |
761 |
|
|
logical :: doesit |
762 |
|
|
doesit = FF_uses_dipoles .or. FF_uses_sticky .or. & |
763 |
|
|
FF_uses_GB .or. FF_uses_RF |
764 |
|
|
end function FF_UsesDirectionalAtoms |
765 |
|
|
|
766 |
|
|
function FF_RequiresPrepairCalc() result(doesit) |
767 |
|
|
logical :: doesit |
768 |
|
|
doesit = FF_uses_EAM |
769 |
|
|
end function FF_RequiresPrepairCalc |
770 |
|
|
|
771 |
|
|
function FF_RequiresPostpairCalc() result(doesit) |
772 |
|
|
logical :: doesit |
773 |
|
|
doesit = FF_uses_RF |
774 |
|
|
end function FF_RequiresPostpairCalc |
775 |
|
|
|
776 |
|
|
end module do_Forces |