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 |
650 |
!! @version $Id: do_Forces.F90,v 1.25 2003-07-24 19:57:35 chuckv Exp $, $Date: 2003-07-24 19:57:35 $, $Name: not supported by cvs2svn $, $Revision: 1.25 $ |
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 |
mmeineke |
626 |
use vector_class |
21 |
chuckv |
650 |
use eam |
22 |
mmeineke |
377 |
#ifdef IS_MPI |
23 |
|
|
use mpiSimulation |
24 |
|
|
#endif |
25 |
|
|
|
26 |
|
|
implicit none |
27 |
|
|
PRIVATE |
28 |
|
|
|
29 |
|
|
#define __FORTRAN90 |
30 |
|
|
#include "fForceField.h" |
31 |
|
|
|
32 |
mmeineke |
626 |
logical, save :: do_forces_initialized = .false., haveRlist = .false. |
33 |
|
|
logical, save :: havePolicies = .false. |
34 |
mmeineke |
377 |
logical, save :: FF_uses_LJ |
35 |
|
|
logical, save :: FF_uses_sticky |
36 |
|
|
logical, save :: FF_uses_dipoles |
37 |
|
|
logical, save :: FF_uses_RF |
38 |
|
|
logical, save :: FF_uses_GB |
39 |
|
|
logical, save :: FF_uses_EAM |
40 |
|
|
|
41 |
mmeineke |
626 |
real(kind=dp), save :: rlist, rlistsq |
42 |
|
|
|
43 |
mmeineke |
377 |
public :: init_FF |
44 |
|
|
public :: do_force_loop |
45 |
mmeineke |
626 |
public :: setRlistDF |
46 |
mmeineke |
377 |
|
47 |
|
|
contains |
48 |
|
|
|
49 |
mmeineke |
626 |
subroutine setRlistDF( this_rlist ) |
50 |
|
|
|
51 |
|
|
real(kind=dp) :: this_rlist |
52 |
|
|
|
53 |
|
|
rlist = this_rlist |
54 |
|
|
rlistsq = rlist * rlist |
55 |
|
|
|
56 |
|
|
haveRlist = .true. |
57 |
|
|
if( havePolicies ) do_forces_initialized = .true. |
58 |
|
|
|
59 |
|
|
end subroutine setRlistDF |
60 |
|
|
|
61 |
mmeineke |
377 |
subroutine init_FF(LJMIXPOLICY, use_RF_c, thisStat) |
62 |
|
|
|
63 |
|
|
integer, intent(in) :: LJMIXPOLICY |
64 |
|
|
logical, intent(in) :: use_RF_c |
65 |
|
|
|
66 |
|
|
integer, intent(out) :: thisStat |
67 |
|
|
integer :: my_status, nMatches |
68 |
|
|
integer, pointer :: MatchList(:) => null() |
69 |
|
|
real(kind=dp) :: rcut, rrf, rt, dielect |
70 |
|
|
|
71 |
|
|
!! assume things are copacetic, unless they aren't |
72 |
|
|
thisStat = 0 |
73 |
|
|
|
74 |
|
|
!! Fortran's version of a cast: |
75 |
|
|
FF_uses_RF = use_RF_c |
76 |
|
|
|
77 |
|
|
!! init_FF is called *after* all of the atom types have been |
78 |
|
|
!! defined in atype_module using the new_atype subroutine. |
79 |
|
|
!! |
80 |
|
|
!! this will scan through the known atypes and figure out what |
81 |
|
|
!! interactions are used by the force field. |
82 |
|
|
|
83 |
|
|
FF_uses_LJ = .false. |
84 |
|
|
FF_uses_sticky = .false. |
85 |
|
|
FF_uses_dipoles = .false. |
86 |
|
|
FF_uses_GB = .false. |
87 |
|
|
FF_uses_EAM = .false. |
88 |
|
|
|
89 |
|
|
call getMatchingElementList(atypes, "is_LJ", .true., nMatches, MatchList) |
90 |
|
|
if (nMatches .gt. 0) FF_uses_LJ = .true. |
91 |
|
|
|
92 |
|
|
call getMatchingElementList(atypes, "is_DP", .true., nMatches, MatchList) |
93 |
|
|
if (nMatches .gt. 0) FF_uses_dipoles = .true. |
94 |
|
|
|
95 |
|
|
call getMatchingElementList(atypes, "is_Sticky", .true., nMatches, & |
96 |
|
|
MatchList) |
97 |
|
|
if (nMatches .gt. 0) FF_uses_Sticky = .true. |
98 |
|
|
|
99 |
|
|
call getMatchingElementList(atypes, "is_GB", .true., nMatches, MatchList) |
100 |
|
|
if (nMatches .gt. 0) FF_uses_GB = .true. |
101 |
|
|
|
102 |
|
|
call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList) |
103 |
|
|
if (nMatches .gt. 0) FF_uses_EAM = .true. |
104 |
|
|
|
105 |
|
|
!! check to make sure the FF_uses_RF setting makes sense |
106 |
|
|
|
107 |
|
|
if (FF_uses_dipoles) then |
108 |
|
|
if (FF_uses_RF) then |
109 |
|
|
dielect = getDielect() |
110 |
mmeineke |
626 |
call initialize_rf(dielect) |
111 |
mmeineke |
377 |
endif |
112 |
|
|
else |
113 |
|
|
if (FF_uses_RF) then |
114 |
|
|
write(default_error,*) 'Using Reaction Field with no dipoles? Huh?' |
115 |
|
|
thisStat = -1 |
116 |
|
|
return |
117 |
|
|
endif |
118 |
mmeineke |
626 |
endif |
119 |
mmeineke |
377 |
|
120 |
|
|
if (FF_uses_LJ) then |
121 |
|
|
|
122 |
|
|
select case (LJMIXPOLICY) |
123 |
|
|
case (LB_MIXING_RULE) |
124 |
mmeineke |
626 |
call init_lj_FF(LB_MIXING_RULE, my_status) |
125 |
mmeineke |
377 |
case (EXPLICIT_MIXING_RULE) |
126 |
mmeineke |
626 |
call init_lj_FF(EXPLICIT_MIXING_RULE, my_status) |
127 |
mmeineke |
377 |
case default |
128 |
|
|
write(default_error,*) 'unknown LJ Mixing Policy!' |
129 |
|
|
thisStat = -1 |
130 |
|
|
return |
131 |
|
|
end select |
132 |
|
|
if (my_status /= 0) then |
133 |
|
|
thisStat = -1 |
134 |
|
|
return |
135 |
|
|
end if |
136 |
|
|
endif |
137 |
|
|
|
138 |
|
|
if (FF_uses_sticky) then |
139 |
|
|
call check_sticky_FF(my_status) |
140 |
|
|
if (my_status /= 0) then |
141 |
|
|
thisStat = -1 |
142 |
|
|
return |
143 |
|
|
end if |
144 |
|
|
endif |
145 |
|
|
|
146 |
|
|
if (FF_uses_GB) then |
147 |
|
|
call check_gb_pair_FF(my_status) |
148 |
|
|
if (my_status .ne. 0) then |
149 |
|
|
thisStat = -1 |
150 |
|
|
return |
151 |
|
|
endif |
152 |
|
|
endif |
153 |
|
|
|
154 |
|
|
if (FF_uses_GB .and. FF_uses_LJ) then |
155 |
|
|
endif |
156 |
chuckv |
480 |
if (.not. do_forces_initialized) then |
157 |
|
|
!! Create neighbor lists |
158 |
|
|
call expandNeighborList(getNlocal(), my_status) |
159 |
|
|
if (my_Status /= 0) then |
160 |
|
|
write(default_error,*) "SimSetup: ExpandNeighborList returned error." |
161 |
|
|
thisStat = -1 |
162 |
|
|
return |
163 |
|
|
endif |
164 |
|
|
endif |
165 |
mmeineke |
377 |
|
166 |
mmeineke |
626 |
havePolicies = .true. |
167 |
|
|
if( haveRlist ) do_forces_initialized = .true. |
168 |
|
|
|
169 |
mmeineke |
377 |
end subroutine init_FF |
170 |
|
|
|
171 |
|
|
|
172 |
|
|
!! Does force loop over i,j pairs. Calls do_pair to calculates forces. |
173 |
|
|
!-------------------------------------------------------------> |
174 |
|
|
subroutine do_force_loop(q, A, u_l, f, t, tau, pot, do_pot_c, do_stress_c, & |
175 |
|
|
error) |
176 |
|
|
!! Position array provided by C, dimensioned by getNlocal |
177 |
|
|
real ( kind = dp ), dimension(3,getNlocal()) :: q |
178 |
|
|
!! Rotation Matrix for each long range particle in simulation. |
179 |
|
|
real( kind = dp), dimension(9,getNlocal()) :: A |
180 |
|
|
!! Unit vectors for dipoles (lab frame) |
181 |
|
|
real( kind = dp ), dimension(3,getNlocal()) :: u_l |
182 |
|
|
!! Force array provided by C, dimensioned by getNlocal |
183 |
|
|
real ( kind = dp ), dimension(3,getNlocal()) :: f |
184 |
|
|
!! Torsion array provided by C, dimensioned by getNlocal |
185 |
|
|
real( kind = dp ), dimension(3,getNlocal()) :: t |
186 |
|
|
!! Stress Tensor |
187 |
|
|
real( kind = dp), dimension(9) :: tau |
188 |
|
|
real ( kind = dp ) :: pot |
189 |
|
|
logical ( kind = 2) :: do_pot_c, do_stress_c |
190 |
|
|
logical :: do_pot |
191 |
|
|
logical :: do_stress |
192 |
chuckv |
439 |
#ifdef IS_MPI |
193 |
chuckv |
441 |
real( kind = DP ) :: pot_local |
194 |
mmeineke |
377 |
integer :: nrow |
195 |
|
|
integer :: ncol |
196 |
|
|
#endif |
197 |
|
|
integer :: nlocal |
198 |
|
|
integer :: natoms |
199 |
|
|
logical :: update_nlist |
200 |
|
|
integer :: i, j, jbeg, jend, jnab |
201 |
|
|
integer :: nlist |
202 |
mmeineke |
626 |
real( kind = DP ) :: rijsq |
203 |
mmeineke |
377 |
real(kind=dp),dimension(3) :: d |
204 |
|
|
real(kind=dp) :: rfpot, mu_i, virial |
205 |
|
|
integer :: me_i |
206 |
|
|
logical :: is_dp_i |
207 |
|
|
integer :: neighborListSize |
208 |
|
|
integer :: listerror, error |
209 |
|
|
integer :: localError |
210 |
mmeineke |
626 |
|
211 |
|
|
real(kind=dp) :: listSkin = 1.0 |
212 |
mmeineke |
619 |
|
213 |
mmeineke |
377 |
|
214 |
|
|
!! initialize local variables |
215 |
|
|
|
216 |
|
|
#ifdef IS_MPI |
217 |
chuckv |
441 |
pot_local = 0.0_dp |
218 |
mmeineke |
377 |
nlocal = getNlocal() |
219 |
|
|
nrow = getNrow(plan_row) |
220 |
|
|
ncol = getNcol(plan_col) |
221 |
|
|
#else |
222 |
|
|
nlocal = getNlocal() |
223 |
|
|
natoms = nlocal |
224 |
|
|
#endif |
225 |
chuckv |
441 |
|
226 |
mmeineke |
377 |
call check_initialization(localError) |
227 |
|
|
if ( localError .ne. 0 ) then |
228 |
|
|
error = -1 |
229 |
|
|
return |
230 |
|
|
end if |
231 |
|
|
call zero_work_arrays() |
232 |
|
|
|
233 |
|
|
do_pot = do_pot_c |
234 |
|
|
do_stress = do_stress_c |
235 |
|
|
|
236 |
|
|
! Gather all information needed by all force loops: |
237 |
|
|
|
238 |
|
|
#ifdef IS_MPI |
239 |
|
|
|
240 |
|
|
call gather(q,q_Row,plan_row3d) |
241 |
|
|
call gather(q,q_Col,plan_col3d) |
242 |
|
|
|
243 |
|
|
if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then |
244 |
|
|
call gather(u_l,u_l_Row,plan_row3d) |
245 |
|
|
call gather(u_l,u_l_Col,plan_col3d) |
246 |
|
|
|
247 |
|
|
call gather(A,A_Row,plan_row_rotation) |
248 |
|
|
call gather(A,A_Col,plan_col_rotation) |
249 |
|
|
endif |
250 |
|
|
|
251 |
|
|
#endif |
252 |
|
|
|
253 |
|
|
if (FF_RequiresPrepairCalc() .and. SimRequiresPrepairCalc()) then |
254 |
|
|
!! See if we need to update neighbor lists |
255 |
mmeineke |
626 |
call checkNeighborList(nlocal, q, listSkin, update_nlist) |
256 |
mmeineke |
377 |
!! if_mpi_gather_stuff_for_prepair |
257 |
|
|
!! do_prepair_loop_if_needed |
258 |
|
|
!! if_mpi_scatter_stuff_from_prepair |
259 |
|
|
!! if_mpi_gather_stuff_from_prepair_to_main_loop |
260 |
chuckv |
648 |
|
261 |
|
|
!--------------------PREFORCE LOOP----------->>>>>>>>>>>>>>>>>>>>>>>>>>> |
262 |
|
|
#ifdef IS_MPI |
263 |
|
|
|
264 |
|
|
if (update_nlist) then |
265 |
|
|
|
266 |
|
|
!! save current configuration, construct neighbor list, |
267 |
|
|
!! and calculate forces |
268 |
|
|
call saveNeighborList(nlocal, q) |
269 |
|
|
|
270 |
|
|
neighborListSize = size(list) |
271 |
|
|
nlist = 0 |
272 |
|
|
|
273 |
|
|
do i = 1, nrow |
274 |
|
|
point(i) = nlist + 1 |
275 |
|
|
|
276 |
|
|
prepair_inner: do j = 1, ncol |
277 |
|
|
|
278 |
|
|
if (skipThisPair(i,j)) cycle prepair_inner |
279 |
|
|
|
280 |
|
|
call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq) |
281 |
|
|
|
282 |
|
|
if (rijsq < rlistsq) then |
283 |
|
|
|
284 |
|
|
nlist = nlist + 1 |
285 |
|
|
|
286 |
|
|
if (nlist > neighborListSize) then |
287 |
|
|
call expandNeighborList(nlocal, listerror) |
288 |
|
|
if (listerror /= 0) then |
289 |
|
|
error = -1 |
290 |
|
|
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
291 |
|
|
return |
292 |
|
|
end if |
293 |
|
|
neighborListSize = size(list) |
294 |
|
|
endif |
295 |
|
|
|
296 |
|
|
list(nlist) = j |
297 |
|
|
call do_prepair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot_local) |
298 |
|
|
endif |
299 |
|
|
enddo prepair_inner |
300 |
|
|
enddo |
301 |
|
|
|
302 |
|
|
point(nrow + 1) = nlist + 1 |
303 |
|
|
|
304 |
|
|
else !! (of update_check) |
305 |
|
|
|
306 |
|
|
! use the list to find the neighbors |
307 |
|
|
do i = 1, nrow |
308 |
|
|
JBEG = POINT(i) |
309 |
|
|
JEND = POINT(i+1) - 1 |
310 |
|
|
! check thiat molecule i has neighbors |
311 |
|
|
if (jbeg .le. jend) then |
312 |
|
|
|
313 |
|
|
do jnab = jbeg, jend |
314 |
|
|
j = list(jnab) |
315 |
|
|
|
316 |
|
|
call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq) |
317 |
|
|
call do_prepair(i, j, rijsq, d, do_pot, do_stress, & |
318 |
|
|
u_l, A, f, t, pot_local) |
319 |
|
|
|
320 |
|
|
enddo |
321 |
|
|
endif |
322 |
|
|
enddo |
323 |
|
|
endif |
324 |
|
|
|
325 |
|
|
#else |
326 |
|
|
|
327 |
|
|
if (update_nlist) then |
328 |
|
|
|
329 |
|
|
! save current configuration, contruct neighbor list, |
330 |
|
|
! and calculate forces |
331 |
|
|
call saveNeighborList(natoms, q) |
332 |
|
|
|
333 |
|
|
neighborListSize = size(list) |
334 |
|
|
|
335 |
|
|
nlist = 0 |
336 |
|
|
|
337 |
|
|
do i = 1, natoms-1 |
338 |
|
|
point(i) = nlist + 1 |
339 |
|
|
|
340 |
|
|
prepair_inner: do j = i+1, natoms |
341 |
|
|
|
342 |
|
|
if (skipThisPair(i,j)) cycle prepair_inner |
343 |
|
|
|
344 |
|
|
call get_interatomic_vector(q(:,i), q(:,j), d, rijsq) |
345 |
|
|
|
346 |
|
|
|
347 |
|
|
if (rijsq < rlistsq) then |
348 |
|
|
|
349 |
|
|
nlist = nlist + 1 |
350 |
|
|
|
351 |
|
|
if (nlist > neighborListSize) then |
352 |
|
|
call expandNeighborList(natoms, listerror) |
353 |
|
|
if (listerror /= 0) then |
354 |
|
|
error = -1 |
355 |
|
|
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
356 |
|
|
return |
357 |
|
|
end if |
358 |
|
|
neighborListSize = size(list) |
359 |
|
|
endif |
360 |
|
|
|
361 |
|
|
list(nlist) = j |
362 |
|
|
|
363 |
|
|
call do_prepair(i, j, rijsq, d, do_pot, do_stress, & |
364 |
|
|
u_l, A, f, t, pot) |
365 |
|
|
|
366 |
|
|
endif |
367 |
|
|
enddo prepair_inner |
368 |
|
|
enddo |
369 |
|
|
|
370 |
|
|
point(natoms) = nlist + 1 |
371 |
|
|
|
372 |
|
|
else !! (update) |
373 |
|
|
|
374 |
|
|
! use the list to find the neighbors |
375 |
|
|
do i = 1, natoms-1 |
376 |
|
|
JBEG = POINT(i) |
377 |
|
|
JEND = POINT(i+1) - 1 |
378 |
|
|
! check thiat molecule i has neighbors |
379 |
|
|
if (jbeg .le. jend) then |
380 |
|
|
|
381 |
|
|
do jnab = jbeg, jend |
382 |
|
|
j = list(jnab) |
383 |
|
|
|
384 |
|
|
call get_interatomic_vector(q(:,i), q(:,j), d, rijsq) |
385 |
|
|
call do_prepair(i, j, rijsq, d, do_pot, do_stress, & |
386 |
|
|
u_l, A, f, t, pot) |
387 |
|
|
|
388 |
|
|
enddo |
389 |
|
|
endif |
390 |
|
|
enddo |
391 |
|
|
endif |
392 |
|
|
#endif |
393 |
|
|
!! Do rest of preforce calculations |
394 |
|
|
call do_preforce(nlocal,pot) |
395 |
mmeineke |
377 |
else |
396 |
chuckv |
648 |
!! See if we need to update neighbor lists for non pre-pair |
397 |
mmeineke |
626 |
call checkNeighborList(nlocal, q, listSkin, update_nlist) |
398 |
mmeineke |
377 |
endif |
399 |
chuckv |
648 |
|
400 |
|
|
|
401 |
|
|
|
402 |
|
|
|
403 |
|
|
|
404 |
|
|
!---------------------------------MAIN Pair LOOP->>>>>>>>>>>>>>>>>>>>>>>>>>>> |
405 |
|
|
|
406 |
|
|
|
407 |
|
|
|
408 |
|
|
|
409 |
|
|
|
410 |
mmeineke |
377 |
#ifdef IS_MPI |
411 |
|
|
|
412 |
|
|
if (update_nlist) then |
413 |
|
|
|
414 |
|
|
!! save current configuration, construct neighbor list, |
415 |
|
|
!! and calculate forces |
416 |
mmeineke |
459 |
call saveNeighborList(nlocal, q) |
417 |
mmeineke |
377 |
|
418 |
|
|
neighborListSize = size(list) |
419 |
|
|
nlist = 0 |
420 |
|
|
|
421 |
|
|
do i = 1, nrow |
422 |
|
|
point(i) = nlist + 1 |
423 |
|
|
|
424 |
|
|
inner: do j = 1, ncol |
425 |
|
|
|
426 |
|
|
if (skipThisPair(i,j)) cycle inner |
427 |
|
|
|
428 |
|
|
call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq) |
429 |
|
|
|
430 |
mmeineke |
626 |
if (rijsq < rlistsq) then |
431 |
mmeineke |
377 |
|
432 |
|
|
nlist = nlist + 1 |
433 |
|
|
|
434 |
|
|
if (nlist > neighborListSize) then |
435 |
|
|
call expandNeighborList(nlocal, listerror) |
436 |
|
|
if (listerror /= 0) then |
437 |
|
|
error = -1 |
438 |
|
|
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
439 |
|
|
return |
440 |
|
|
end if |
441 |
|
|
neighborListSize = size(list) |
442 |
|
|
endif |
443 |
|
|
|
444 |
|
|
list(nlist) = j |
445 |
|
|
|
446 |
mmeineke |
626 |
call do_pair(i, j, rijsq, d, do_pot, do_stress, & |
447 |
|
|
u_l, A, f, t, pot_local) |
448 |
|
|
|
449 |
mmeineke |
377 |
endif |
450 |
|
|
enddo inner |
451 |
|
|
enddo |
452 |
|
|
|
453 |
|
|
point(nrow + 1) = nlist + 1 |
454 |
|
|
|
455 |
|
|
else !! (of update_check) |
456 |
|
|
|
457 |
|
|
! use the list to find the neighbors |
458 |
|
|
do i = 1, nrow |
459 |
|
|
JBEG = POINT(i) |
460 |
|
|
JEND = POINT(i+1) - 1 |
461 |
|
|
! check thiat molecule i has neighbors |
462 |
|
|
if (jbeg .le. jend) then |
463 |
|
|
|
464 |
|
|
do jnab = jbeg, jend |
465 |
|
|
j = list(jnab) |
466 |
|
|
|
467 |
|
|
call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq) |
468 |
|
|
call do_pair(i, j, rijsq, d, do_pot, do_stress, & |
469 |
chuckv |
441 |
u_l, A, f, t, pot_local) |
470 |
mmeineke |
377 |
|
471 |
|
|
enddo |
472 |
|
|
endif |
473 |
|
|
enddo |
474 |
|
|
endif |
475 |
|
|
|
476 |
|
|
#else |
477 |
|
|
|
478 |
|
|
if (update_nlist) then |
479 |
|
|
|
480 |
|
|
! save current configuration, contruct neighbor list, |
481 |
|
|
! and calculate forces |
482 |
mmeineke |
459 |
call saveNeighborList(natoms, q) |
483 |
mmeineke |
377 |
|
484 |
|
|
neighborListSize = size(list) |
485 |
|
|
|
486 |
|
|
nlist = 0 |
487 |
|
|
|
488 |
|
|
do i = 1, natoms-1 |
489 |
|
|
point(i) = nlist + 1 |
490 |
|
|
|
491 |
|
|
inner: do j = i+1, natoms |
492 |
|
|
|
493 |
chuckv |
388 |
if (skipThisPair(i,j)) cycle inner |
494 |
|
|
|
495 |
mmeineke |
377 |
call get_interatomic_vector(q(:,i), q(:,j), d, rijsq) |
496 |
|
|
|
497 |
chuckv |
388 |
|
498 |
mmeineke |
626 |
if (rijsq < rlistsq) then |
499 |
mmeineke |
377 |
|
500 |
|
|
nlist = nlist + 1 |
501 |
|
|
|
502 |
|
|
if (nlist > neighborListSize) then |
503 |
|
|
call expandNeighborList(natoms, listerror) |
504 |
|
|
if (listerror /= 0) then |
505 |
|
|
error = -1 |
506 |
|
|
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
507 |
|
|
return |
508 |
|
|
end if |
509 |
|
|
neighborListSize = size(list) |
510 |
|
|
endif |
511 |
|
|
|
512 |
|
|
list(nlist) = j |
513 |
|
|
|
514 |
mmeineke |
626 |
call do_pair(i, j, rijsq, d, do_pot, do_stress, & |
515 |
chuckv |
441 |
u_l, A, f, t, pot) |
516 |
mmeineke |
626 |
|
517 |
mmeineke |
377 |
endif |
518 |
|
|
enddo inner |
519 |
|
|
enddo |
520 |
|
|
|
521 |
|
|
point(natoms) = nlist + 1 |
522 |
|
|
|
523 |
|
|
else !! (update) |
524 |
|
|
|
525 |
|
|
! use the list to find the neighbors |
526 |
|
|
do i = 1, natoms-1 |
527 |
|
|
JBEG = POINT(i) |
528 |
|
|
JEND = POINT(i+1) - 1 |
529 |
|
|
! check thiat molecule i has neighbors |
530 |
|
|
if (jbeg .le. jend) then |
531 |
|
|
|
532 |
|
|
do jnab = jbeg, jend |
533 |
|
|
j = list(jnab) |
534 |
|
|
|
535 |
|
|
call get_interatomic_vector(q(:,i), q(:,j), d, rijsq) |
536 |
|
|
call do_pair(i, j, rijsq, d, do_pot, do_stress, & |
537 |
chuckv |
441 |
u_l, A, f, t, pot) |
538 |
mmeineke |
377 |
|
539 |
|
|
enddo |
540 |
|
|
endif |
541 |
|
|
enddo |
542 |
|
|
endif |
543 |
|
|
|
544 |
|
|
#endif |
545 |
|
|
|
546 |
|
|
! phew, done with main loop. |
547 |
|
|
|
548 |
|
|
#ifdef IS_MPI |
549 |
|
|
!!distribute forces |
550 |
chuckv |
438 |
|
551 |
|
|
f_temp = 0.0_dp |
552 |
|
|
call scatter(f_Row,f_temp,plan_row3d) |
553 |
|
|
do i = 1,nlocal |
554 |
|
|
f(1:3,i) = f(1:3,i) + f_temp(1:3,i) |
555 |
|
|
end do |
556 |
|
|
|
557 |
|
|
f_temp = 0.0_dp |
558 |
mmeineke |
377 |
call scatter(f_Col,f_temp,plan_col3d) |
559 |
|
|
do i = 1,nlocal |
560 |
|
|
f(1:3,i) = f(1:3,i) + f_temp(1:3,i) |
561 |
|
|
end do |
562 |
|
|
|
563 |
|
|
if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then |
564 |
chuckv |
438 |
t_temp = 0.0_dp |
565 |
|
|
call scatter(t_Row,t_temp,plan_row3d) |
566 |
|
|
do i = 1,nlocal |
567 |
|
|
t(1:3,i) = t(1:3,i) + t_temp(1:3,i) |
568 |
|
|
end do |
569 |
|
|
t_temp = 0.0_dp |
570 |
mmeineke |
377 |
call scatter(t_Col,t_temp,plan_col3d) |
571 |
|
|
|
572 |
|
|
do i = 1,nlocal |
573 |
|
|
t(1:3,i) = t(1:3,i) + t_temp(1:3,i) |
574 |
|
|
end do |
575 |
|
|
endif |
576 |
|
|
|
577 |
|
|
if (do_pot) then |
578 |
|
|
! scatter/gather pot_row into the members of my column |
579 |
|
|
call scatter(pot_Row, pot_Temp, plan_row) |
580 |
chuckv |
439 |
|
581 |
mmeineke |
377 |
! scatter/gather pot_local into all other procs |
582 |
|
|
! add resultant to get total pot |
583 |
|
|
do i = 1, nlocal |
584 |
|
|
pot_local = pot_local + pot_Temp(i) |
585 |
|
|
enddo |
586 |
chuckv |
439 |
|
587 |
|
|
pot_Temp = 0.0_DP |
588 |
mmeineke |
377 |
|
589 |
|
|
call scatter(pot_Col, pot_Temp, plan_col) |
590 |
|
|
do i = 1, nlocal |
591 |
|
|
pot_local = pot_local + pot_Temp(i) |
592 |
|
|
enddo |
593 |
chuckv |
439 |
|
594 |
mmeineke |
377 |
endif |
595 |
|
|
#endif |
596 |
|
|
|
597 |
|
|
if (FF_RequiresPostpairCalc() .and. SimRequiresPostpairCalc()) then |
598 |
|
|
|
599 |
|
|
if (FF_uses_RF .and. SimUsesRF()) then |
600 |
|
|
|
601 |
|
|
#ifdef IS_MPI |
602 |
|
|
call scatter(rf_Row,rf,plan_row3d) |
603 |
|
|
call scatter(rf_Col,rf_Temp,plan_col3d) |
604 |
|
|
do i = 1,nlocal |
605 |
|
|
rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i) |
606 |
|
|
end do |
607 |
|
|
#endif |
608 |
|
|
|
609 |
|
|
do i = 1, getNlocal() |
610 |
|
|
|
611 |
|
|
rfpot = 0.0_DP |
612 |
|
|
#ifdef IS_MPI |
613 |
|
|
me_i = atid_row(i) |
614 |
|
|
#else |
615 |
|
|
me_i = atid(i) |
616 |
|
|
#endif |
617 |
|
|
call getElementProperty(atypes, me_i, "is_DP", is_DP_i) |
618 |
|
|
if ( is_DP_i ) then |
619 |
|
|
call getElementProperty(atypes, me_i, "dipole_moment", mu_i) |
620 |
|
|
!! The reaction field needs to include a self contribution |
621 |
|
|
!! to the field: |
622 |
|
|
call accumulate_self_rf(i, mu_i, u_l) |
623 |
|
|
!! Get the reaction field contribution to the |
624 |
|
|
!! potential and torques: |
625 |
|
|
call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot) |
626 |
|
|
#ifdef IS_MPI |
627 |
|
|
pot_local = pot_local + rfpot |
628 |
|
|
#else |
629 |
|
|
pot = pot + rfpot |
630 |
|
|
|
631 |
|
|
#endif |
632 |
|
|
endif |
633 |
|
|
enddo |
634 |
|
|
endif |
635 |
|
|
endif |
636 |
|
|
|
637 |
|
|
|
638 |
|
|
#ifdef IS_MPI |
639 |
|
|
|
640 |
|
|
if (do_pot) then |
641 |
chuckv |
441 |
pot = pot + pot_local |
642 |
mmeineke |
377 |
!! we assume the c code will do the allreduce to get the total potential |
643 |
|
|
!! we could do it right here if we needed to... |
644 |
|
|
endif |
645 |
|
|
|
646 |
|
|
if (do_stress) then |
647 |
gezelter |
490 |
call mpi_allreduce(tau_Temp, tau, 9,mpi_double_precision,mpi_sum, & |
648 |
mmeineke |
377 |
mpi_comm_world,mpi_err) |
649 |
chuckv |
470 |
call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, & |
650 |
mmeineke |
377 |
mpi_comm_world,mpi_err) |
651 |
|
|
endif |
652 |
|
|
|
653 |
|
|
#else |
654 |
|
|
|
655 |
|
|
if (do_stress) then |
656 |
|
|
tau = tau_Temp |
657 |
|
|
virial = virial_Temp |
658 |
|
|
endif |
659 |
|
|
|
660 |
|
|
#endif |
661 |
|
|
|
662 |
|
|
end subroutine do_force_loop |
663 |
|
|
|
664 |
chuckv |
441 |
subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot) |
665 |
mmeineke |
377 |
|
666 |
|
|
real( kind = dp ) :: pot |
667 |
chuckv |
460 |
real( kind = dp ), dimension(3,getNlocal()) :: u_l |
668 |
|
|
real (kind=dp), dimension(9,getNlocal()) :: A |
669 |
|
|
real (kind=dp), dimension(3,getNlocal()) :: f |
670 |
|
|
real (kind=dp), dimension(3,getNlocal()) :: t |
671 |
mmeineke |
377 |
|
672 |
|
|
logical, intent(inout) :: do_pot, do_stress |
673 |
|
|
integer, intent(in) :: i, j |
674 |
|
|
real ( kind = dp ), intent(inout) :: rijsq |
675 |
|
|
real ( kind = dp ) :: r |
676 |
|
|
real ( kind = dp ), intent(inout) :: d(3) |
677 |
|
|
logical :: is_LJ_i, is_LJ_j |
678 |
|
|
logical :: is_DP_i, is_DP_j |
679 |
|
|
logical :: is_GB_i, is_GB_j |
680 |
chuckv |
648 |
logical :: is_EAM_i,is_EAM_j |
681 |
mmeineke |
377 |
logical :: is_Sticky_i, is_Sticky_j |
682 |
|
|
integer :: me_i, me_j |
683 |
|
|
|
684 |
|
|
r = sqrt(rijsq) |
685 |
|
|
|
686 |
|
|
#ifdef IS_MPI |
687 |
gezelter |
490 |
if (tagRow(i) .eq. tagColumn(j)) then |
688 |
|
|
write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j) |
689 |
|
|
endif |
690 |
mmeineke |
377 |
|
691 |
|
|
me_i = atid_row(i) |
692 |
|
|
me_j = atid_col(j) |
693 |
|
|
|
694 |
|
|
#else |
695 |
|
|
|
696 |
|
|
me_i = atid(i) |
697 |
|
|
me_j = atid(j) |
698 |
|
|
|
699 |
|
|
#endif |
700 |
|
|
|
701 |
|
|
if (FF_uses_LJ .and. SimUsesLJ()) then |
702 |
|
|
call getElementProperty(atypes, me_i, "is_LJ", is_LJ_i) |
703 |
|
|
call getElementProperty(atypes, me_j, "is_LJ", is_LJ_j) |
704 |
|
|
|
705 |
|
|
if ( is_LJ_i .and. is_LJ_j ) & |
706 |
|
|
call do_lj_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress) |
707 |
|
|
endif |
708 |
|
|
|
709 |
|
|
if (FF_uses_dipoles .and. SimUsesDipoles()) then |
710 |
|
|
call getElementProperty(atypes, me_i, "is_DP", is_DP_i) |
711 |
|
|
call getElementProperty(atypes, me_j, "is_DP", is_DP_j) |
712 |
|
|
|
713 |
|
|
if ( is_DP_i .and. is_DP_j ) then |
714 |
|
|
|
715 |
gezelter |
462 |
call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, & |
716 |
mmeineke |
377 |
do_pot, do_stress) |
717 |
|
|
if (FF_uses_RF .and. SimUsesRF()) then |
718 |
|
|
call accumulate_rf(i, j, r, u_l) |
719 |
|
|
call rf_correct_forces(i, j, d, r, u_l, f, do_stress) |
720 |
|
|
endif |
721 |
|
|
|
722 |
|
|
endif |
723 |
|
|
endif |
724 |
|
|
|
725 |
|
|
if (FF_uses_Sticky .and. SimUsesSticky()) then |
726 |
|
|
|
727 |
|
|
call getElementProperty(atypes, me_i, "is_Sticky", is_Sticky_i) |
728 |
|
|
call getElementProperty(atypes, me_j, "is_Sticky", is_Sticky_j) |
729 |
chuckv |
388 |
|
730 |
mmeineke |
377 |
if ( is_Sticky_i .and. is_Sticky_j ) then |
731 |
|
|
call do_sticky_pair(i, j, d, r, rijsq, A, pot, f, t, & |
732 |
|
|
do_pot, do_stress) |
733 |
|
|
endif |
734 |
|
|
endif |
735 |
|
|
|
736 |
|
|
|
737 |
|
|
if (FF_uses_GB .and. SimUsesGB()) then |
738 |
|
|
|
739 |
|
|
call getElementProperty(atypes, me_i, "is_GB", is_GB_i) |
740 |
|
|
call getElementProperty(atypes, me_j, "is_GB", is_GB_j) |
741 |
|
|
|
742 |
|
|
if ( is_GB_i .and. is_GB_j ) then |
743 |
|
|
call do_gb_pair(i, j, d, r, rijsq, u_l, pot, f, t, & |
744 |
|
|
do_pot, do_stress) |
745 |
|
|
endif |
746 |
|
|
endif |
747 |
|
|
|
748 |
mmeineke |
597 |
|
749 |
chuckv |
648 |
|
750 |
|
|
if (FF_uses_EAM .and. SimUsesEAM()) then |
751 |
|
|
call getElementProperty(atypes, me_i, "is_EAM", is_EAM_i) |
752 |
|
|
call getElementProperty(atypes, me_j, "is_EAM", is_EAM_j) |
753 |
|
|
|
754 |
|
|
if ( is_EAM_i .and. is_EAM_j ) & |
755 |
|
|
call do_eam_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress) |
756 |
|
|
endif |
757 |
|
|
|
758 |
mmeineke |
597 |
|
759 |
chuckv |
648 |
|
760 |
|
|
|
761 |
mmeineke |
377 |
end subroutine do_pair |
762 |
|
|
|
763 |
|
|
|
764 |
chuckv |
631 |
|
765 |
chuckv |
648 |
subroutine do_prepair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot) |
766 |
chuckv |
631 |
real( kind = dp ) :: pot |
767 |
|
|
real( kind = dp ), dimension(3,getNlocal()) :: u_l |
768 |
|
|
real (kind=dp), dimension(9,getNlocal()) :: A |
769 |
|
|
real (kind=dp), dimension(3,getNlocal()) :: f |
770 |
|
|
real (kind=dp), dimension(3,getNlocal()) :: t |
771 |
|
|
|
772 |
|
|
logical, intent(inout) :: do_pot, do_stress |
773 |
|
|
integer, intent(in) :: i, j |
774 |
|
|
real ( kind = dp ), intent(inout) :: rijsq |
775 |
|
|
real ( kind = dp ) :: r |
776 |
|
|
real ( kind = dp ), intent(inout) :: d(3) |
777 |
|
|
|
778 |
|
|
logical :: is_EAM_i, is_EAM_j |
779 |
|
|
|
780 |
|
|
integer :: me_i, me_j |
781 |
|
|
|
782 |
|
|
r = sqrt(rijsq) |
783 |
|
|
|
784 |
|
|
#ifdef IS_MPI |
785 |
|
|
if (tagRow(i) .eq. tagColumn(j)) then |
786 |
|
|
write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j) |
787 |
|
|
endif |
788 |
|
|
|
789 |
|
|
me_i = atid_row(i) |
790 |
|
|
me_j = atid_col(j) |
791 |
|
|
|
792 |
|
|
#else |
793 |
|
|
|
794 |
|
|
me_i = atid(i) |
795 |
|
|
me_j = atid(j) |
796 |
|
|
|
797 |
|
|
#endif |
798 |
|
|
|
799 |
|
|
if (FF_uses_EAM .and. SimUsesEAM()) then |
800 |
|
|
call getElementProperty(atypes, me_i, "is_EAM", is_EAM_i) |
801 |
|
|
call getElementProperty(atypes, me_j, "is_EAM", is_EAM_j) |
802 |
|
|
|
803 |
chuckv |
648 |
if ( is_EAM_i .and. is_EAM_j ) & |
804 |
|
|
call calc_EAM_prepair_rho(i, j, d, r, rijsq ) |
805 |
chuckv |
631 |
endif |
806 |
chuckv |
648 |
end subroutine do_prepair |
807 |
|
|
|
808 |
|
|
|
809 |
|
|
|
810 |
|
|
|
811 |
|
|
subroutine do_preforce(nlocal,pot) |
812 |
|
|
integer :: nlocal |
813 |
|
|
real( kind = dp ) :: pot |
814 |
|
|
|
815 |
|
|
if (FF_uses_EAM .and. SimUsesEAM()) then |
816 |
|
|
call calc_EAM_preforce_Frho(nlocal,pot) |
817 |
|
|
endif |
818 |
|
|
|
819 |
|
|
|
820 |
chuckv |
631 |
end subroutine do_preforce |
821 |
|
|
|
822 |
|
|
|
823 |
mmeineke |
377 |
subroutine get_interatomic_vector(q_i, q_j, d, r_sq) |
824 |
|
|
|
825 |
|
|
real (kind = dp), dimension(3) :: q_i |
826 |
|
|
real (kind = dp), dimension(3) :: q_j |
827 |
|
|
real ( kind = dp ), intent(out) :: r_sq |
828 |
gezelter |
571 |
real( kind = dp ) :: d(3), scaled(3) |
829 |
|
|
integer i |
830 |
|
|
|
831 |
chuckv |
482 |
d(1:3) = q_j(1:3) - q_i(1:3) |
832 |
gezelter |
571 |
|
833 |
mmeineke |
377 |
! Wrap back into periodic box if necessary |
834 |
|
|
if ( SimUsesPBC() ) then |
835 |
mmeineke |
393 |
|
836 |
gezelter |
571 |
if( .not.boxIsOrthorhombic ) then |
837 |
|
|
! calc the scaled coordinates. |
838 |
|
|
|
839 |
mmeineke |
572 |
scaled = matmul(HmatInv, d) |
840 |
gezelter |
571 |
|
841 |
|
|
! wrap the scaled coordinates |
842 |
|
|
|
843 |
mmeineke |
572 |
scaled = scaled - anint(scaled) |
844 |
|
|
|
845 |
gezelter |
571 |
|
846 |
|
|
! calc the wrapped real coordinates from the wrapped scaled |
847 |
|
|
! coordinates |
848 |
|
|
|
849 |
mmeineke |
572 |
d = matmul(Hmat,scaled) |
850 |
gezelter |
571 |
|
851 |
|
|
else |
852 |
|
|
! calc the scaled coordinates. |
853 |
|
|
|
854 |
|
|
do i = 1, 3 |
855 |
|
|
scaled(i) = d(i) * HmatInv(i,i) |
856 |
|
|
|
857 |
|
|
! wrap the scaled coordinates |
858 |
|
|
|
859 |
|
|
scaled(i) = scaled(i) - anint(scaled(i)) |
860 |
|
|
|
861 |
|
|
! calc the wrapped real coordinates from the wrapped scaled |
862 |
|
|
! coordinates |
863 |
|
|
|
864 |
|
|
d(i) = scaled(i)*Hmat(i,i) |
865 |
|
|
enddo |
866 |
|
|
endif |
867 |
mmeineke |
393 |
|
868 |
mmeineke |
377 |
endif |
869 |
gezelter |
571 |
|
870 |
mmeineke |
377 |
r_sq = dot_product(d,d) |
871 |
gezelter |
571 |
|
872 |
mmeineke |
377 |
end subroutine get_interatomic_vector |
873 |
gezelter |
571 |
|
874 |
mmeineke |
377 |
subroutine check_initialization(error) |
875 |
|
|
integer, intent(out) :: error |
876 |
|
|
|
877 |
|
|
error = 0 |
878 |
|
|
! Make sure we are properly initialized. |
879 |
|
|
if (.not. do_forces_initialized) then |
880 |
|
|
error = -1 |
881 |
|
|
return |
882 |
|
|
endif |
883 |
|
|
|
884 |
|
|
#ifdef IS_MPI |
885 |
|
|
if (.not. isMPISimSet()) then |
886 |
|
|
write(default_error,*) "ERROR: mpiSimulation has not been initialized!" |
887 |
|
|
error = -1 |
888 |
|
|
return |
889 |
|
|
endif |
890 |
|
|
#endif |
891 |
|
|
|
892 |
|
|
return |
893 |
|
|
end subroutine check_initialization |
894 |
|
|
|
895 |
|
|
|
896 |
|
|
subroutine zero_work_arrays() |
897 |
|
|
|
898 |
|
|
#ifdef IS_MPI |
899 |
|
|
|
900 |
|
|
q_Row = 0.0_dp |
901 |
|
|
q_Col = 0.0_dp |
902 |
|
|
|
903 |
|
|
u_l_Row = 0.0_dp |
904 |
|
|
u_l_Col = 0.0_dp |
905 |
|
|
|
906 |
|
|
A_Row = 0.0_dp |
907 |
|
|
A_Col = 0.0_dp |
908 |
|
|
|
909 |
|
|
f_Row = 0.0_dp |
910 |
|
|
f_Col = 0.0_dp |
911 |
|
|
f_Temp = 0.0_dp |
912 |
|
|
|
913 |
|
|
t_Row = 0.0_dp |
914 |
|
|
t_Col = 0.0_dp |
915 |
|
|
t_Temp = 0.0_dp |
916 |
|
|
|
917 |
|
|
pot_Row = 0.0_dp |
918 |
|
|
pot_Col = 0.0_dp |
919 |
|
|
pot_Temp = 0.0_dp |
920 |
|
|
|
921 |
|
|
rf_Row = 0.0_dp |
922 |
|
|
rf_Col = 0.0_dp |
923 |
|
|
rf_Temp = 0.0_dp |
924 |
|
|
|
925 |
|
|
#endif |
926 |
|
|
|
927 |
|
|
rf = 0.0_dp |
928 |
|
|
tau_Temp = 0.0_dp |
929 |
|
|
virial_Temp = 0.0_dp |
930 |
|
|
end subroutine zero_work_arrays |
931 |
|
|
|
932 |
|
|
function skipThisPair(atom1, atom2) result(skip_it) |
933 |
|
|
integer, intent(in) :: atom1 |
934 |
|
|
integer, intent(in), optional :: atom2 |
935 |
|
|
logical :: skip_it |
936 |
|
|
integer :: unique_id_1, unique_id_2 |
937 |
chuckv |
388 |
integer :: me_i,me_j |
938 |
mmeineke |
377 |
integer :: i |
939 |
|
|
|
940 |
|
|
skip_it = .false. |
941 |
|
|
|
942 |
|
|
!! there are a number of reasons to skip a pair or a particle |
943 |
|
|
!! mostly we do this to exclude atoms who are involved in short |
944 |
|
|
!! range interactions (bonds, bends, torsions), but we also need |
945 |
|
|
!! to exclude some overcounted interactions that result from |
946 |
|
|
!! the parallel decomposition |
947 |
|
|
|
948 |
|
|
#ifdef IS_MPI |
949 |
|
|
!! in MPI, we have to look up the unique IDs for each atom |
950 |
|
|
unique_id_1 = tagRow(atom1) |
951 |
|
|
#else |
952 |
|
|
!! in the normal loop, the atom numbers are unique |
953 |
|
|
unique_id_1 = atom1 |
954 |
|
|
#endif |
955 |
chuckv |
388 |
|
956 |
mmeineke |
377 |
!! We were called with only one atom, so just check the global exclude |
957 |
|
|
!! list for this atom |
958 |
|
|
if (.not. present(atom2)) then |
959 |
|
|
do i = 1, nExcludes_global |
960 |
|
|
if (excludesGlobal(i) == unique_id_1) then |
961 |
|
|
skip_it = .true. |
962 |
|
|
return |
963 |
|
|
end if |
964 |
|
|
end do |
965 |
|
|
return |
966 |
|
|
end if |
967 |
|
|
|
968 |
|
|
#ifdef IS_MPI |
969 |
|
|
unique_id_2 = tagColumn(atom2) |
970 |
|
|
#else |
971 |
|
|
unique_id_2 = atom2 |
972 |
|
|
#endif |
973 |
chuckv |
441 |
|
974 |
mmeineke |
377 |
#ifdef IS_MPI |
975 |
|
|
!! this situation should only arise in MPI simulations |
976 |
|
|
if (unique_id_1 == unique_id_2) then |
977 |
|
|
skip_it = .true. |
978 |
|
|
return |
979 |
|
|
end if |
980 |
|
|
|
981 |
|
|
!! this prevents us from doing the pair on multiple processors |
982 |
|
|
if (unique_id_1 < unique_id_2) then |
983 |
chuckv |
441 |
if (mod(unique_id_1 + unique_id_2,2) == 0) then |
984 |
|
|
skip_it = .true. |
985 |
|
|
return |
986 |
|
|
endif |
987 |
mmeineke |
377 |
else |
988 |
chuckv |
441 |
if (mod(unique_id_1 + unique_id_2,2) == 1) then |
989 |
|
|
skip_it = .true. |
990 |
|
|
return |
991 |
|
|
endif |
992 |
mmeineke |
377 |
endif |
993 |
|
|
#endif |
994 |
chuckv |
441 |
|
995 |
mmeineke |
377 |
!! the rest of these situations can happen in all simulations: |
996 |
|
|
do i = 1, nExcludes_global |
997 |
|
|
if ((excludesGlobal(i) == unique_id_1) .or. & |
998 |
|
|
(excludesGlobal(i) == unique_id_2)) then |
999 |
|
|
skip_it = .true. |
1000 |
|
|
return |
1001 |
|
|
endif |
1002 |
|
|
enddo |
1003 |
chuckv |
441 |
|
1004 |
mmeineke |
377 |
do i = 1, nExcludes_local |
1005 |
|
|
if (excludesLocal(1,i) == unique_id_1) then |
1006 |
|
|
if (excludesLocal(2,i) == unique_id_2) then |
1007 |
|
|
skip_it = .true. |
1008 |
|
|
return |
1009 |
|
|
endif |
1010 |
|
|
else |
1011 |
|
|
if (excludesLocal(1,i) == unique_id_2) then |
1012 |
|
|
if (excludesLocal(2,i) == unique_id_1) then |
1013 |
|
|
skip_it = .true. |
1014 |
|
|
return |
1015 |
|
|
endif |
1016 |
|
|
endif |
1017 |
|
|
endif |
1018 |
|
|
end do |
1019 |
|
|
|
1020 |
|
|
return |
1021 |
|
|
end function skipThisPair |
1022 |
|
|
|
1023 |
|
|
function FF_UsesDirectionalAtoms() result(doesit) |
1024 |
|
|
logical :: doesit |
1025 |
|
|
doesit = FF_uses_dipoles .or. FF_uses_sticky .or. & |
1026 |
|
|
FF_uses_GB .or. FF_uses_RF |
1027 |
|
|
end function FF_UsesDirectionalAtoms |
1028 |
|
|
|
1029 |
|
|
function FF_RequiresPrepairCalc() result(doesit) |
1030 |
|
|
logical :: doesit |
1031 |
|
|
doesit = FF_uses_EAM |
1032 |
|
|
end function FF_RequiresPrepairCalc |
1033 |
|
|
|
1034 |
|
|
function FF_RequiresPostpairCalc() result(doesit) |
1035 |
|
|
logical :: doesit |
1036 |
|
|
doesit = FF_uses_RF |
1037 |
|
|
end function FF_RequiresPostpairCalc |
1038 |
|
|
|
1039 |
|
|
end module do_Forces |