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