ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 619
Committed: Tue Jul 15 22:22:41 2003 UTC (20 years, 11 months ago) by mmeineke
File size: 22678 byte(s)
Log Message:
fixed a long lived bug in do_forces. Rrf was not being used in the neighborlist correctly. rcut was conssistently being set lowere than Rrf causing the dipole cutoff region to be to small. Also led to the removal of the taper region to buffer the dipole cutoff.

File Contents

# Content
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.20 2003-07-15 22:22:41 mmeineke Exp $, $Date: 2003-07-15 22:22:41 $, $Name: not supported by cvs2svn $, $Revision: 1.20 $
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 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 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
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 #ifdef IS_MPI
179 real( kind = DP ) :: pot_local
180 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, rrf, rt, dielect
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
198 !! initialize local variables
199
200 #ifdef IS_MPI
201 pot_local = 0.0_dp
202 nlocal = getNlocal()
203 nrow = getNrow(plan_row)
204 ncol = getNcol(plan_col)
205 #else
206 nlocal = getNlocal()
207 natoms = nlocal
208 #endif
209
210 call getRcut(rcut,rc2=rcutsq)
211 call getRlist(rlist,rlistsq)
212 rt = getRt()
213 rrf = getRrf()
214 dielect = getDielect()
215
216 if( FF_uses_LJ) then
217 call lj_new_rcut( rcut, localError )
218 if ( localError .ne. 0 ) then
219 error = -1
220 return
221 end if
222 end if
223
224
225 if( FF_uses_dipoles ) then
226
227 if( rcut .lt. rrf ) then
228 rcut = rrf
229 rlist = rcut + 1.0_dp
230 rcutsq = rcut * rcut
231 rlistsq = rlist * rlist
232 end if
233
234 call initialize_dipole( rrf, rt )
235 end if
236
237 if( FF_uses_RF )call initialize_rf( rrf, rt, dielect )
238
239
240 call check_initialization(localError)
241 if ( localError .ne. 0 ) then
242 error = -1
243 return
244 end if
245 call zero_work_arrays()
246
247 do_pot = do_pot_c
248 do_stress = do_stress_c
249
250 ! Gather all information needed by all force loops:
251
252 #ifdef IS_MPI
253
254 call gather(q,q_Row,plan_row3d)
255 call gather(q,q_Col,plan_col3d)
256
257 if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then
258 call gather(u_l,u_l_Row,plan_row3d)
259 call gather(u_l,u_l_Col,plan_col3d)
260
261 call gather(A,A_Row,plan_row_rotation)
262 call gather(A,A_Col,plan_col_rotation)
263 endif
264
265 #endif
266
267 if (FF_RequiresPrepairCalc() .and. SimRequiresPrepairCalc()) then
268 !! See if we need to update neighbor lists
269 call checkNeighborList(nlocal, q, rcut, rlist, update_nlist)
270 !! if_mpi_gather_stuff_for_prepair
271 !! do_prepair_loop_if_needed
272 !! if_mpi_scatter_stuff_from_prepair
273 !! if_mpi_gather_stuff_from_prepair_to_main_loop
274 else
275 !! See if we need to update neighbor lists
276 call checkNeighborList(nlocal, q, rcut, rlist, update_nlist)
277 endif
278
279 #ifdef IS_MPI
280
281 if (update_nlist) then
282
283 !! save current configuration, construct neighbor list,
284 !! and calculate forces
285 call saveNeighborList(nlocal, q)
286
287 neighborListSize = size(list)
288 nlist = 0
289
290 do i = 1, nrow
291 point(i) = nlist + 1
292
293 inner: do j = 1, ncol
294
295 if (skipThisPair(i,j)) cycle inner
296
297 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
298
299 if (rijsq < rlistsq) then
300
301 nlist = nlist + 1
302
303 if (nlist > neighborListSize) then
304 call expandNeighborList(nlocal, listerror)
305 if (listerror /= 0) then
306 error = -1
307 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
308 return
309 end if
310 neighborListSize = size(list)
311 endif
312
313 list(nlist) = j
314
315 if (rijsq < rcutsq) then
316 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
317 u_l, A, f, t, pot_local)
318 endif
319 endif
320 enddo inner
321 enddo
322
323 point(nrow + 1) = nlist + 1
324
325 else !! (of update_check)
326
327 ! use the list to find the neighbors
328 do i = 1, nrow
329 JBEG = POINT(i)
330 JEND = POINT(i+1) - 1
331 ! check thiat molecule i has neighbors
332 if (jbeg .le. jend) then
333
334 do jnab = jbeg, jend
335 j = list(jnab)
336
337 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
338 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
339 u_l, A, f, t, pot_local)
340
341 enddo
342 endif
343 enddo
344 endif
345
346 #else
347
348 if (update_nlist) then
349
350 ! save current configuration, contruct neighbor list,
351 ! and calculate forces
352 call saveNeighborList(natoms, q)
353
354 neighborListSize = size(list)
355
356 nlist = 0
357
358 do i = 1, natoms-1
359 point(i) = nlist + 1
360
361 inner: do j = i+1, natoms
362
363 if (skipThisPair(i,j)) cycle inner
364
365 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
366
367
368 if (rijsq < rlistsq) then
369
370 nlist = nlist + 1
371
372 if (nlist > neighborListSize) then
373 call expandNeighborList(natoms, listerror)
374 if (listerror /= 0) then
375 error = -1
376 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
377 return
378 end if
379 neighborListSize = size(list)
380 endif
381
382 list(nlist) = j
383
384 if (rijsq < rcutsq) then
385 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
386 u_l, A, f, t, pot)
387 endif
388 endif
389 enddo inner
390 enddo
391
392 point(natoms) = nlist + 1
393
394 else !! (update)
395
396 ! use the list to find the neighbors
397 do i = 1, natoms-1
398 JBEG = POINT(i)
399 JEND = POINT(i+1) - 1
400 ! check thiat molecule i has neighbors
401 if (jbeg .le. jend) then
402
403 do jnab = jbeg, jend
404 j = list(jnab)
405
406 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
407 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
408 u_l, A, f, t, pot)
409
410 enddo
411 endif
412 enddo
413 endif
414
415 #endif
416
417 ! phew, done with main loop.
418
419 #ifdef IS_MPI
420 !!distribute forces
421
422 f_temp = 0.0_dp
423 call scatter(f_Row,f_temp,plan_row3d)
424 do i = 1,nlocal
425 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
426 end do
427
428 f_temp = 0.0_dp
429 call scatter(f_Col,f_temp,plan_col3d)
430 do i = 1,nlocal
431 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
432 end do
433
434 if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then
435 t_temp = 0.0_dp
436 call scatter(t_Row,t_temp,plan_row3d)
437 do i = 1,nlocal
438 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
439 end do
440 t_temp = 0.0_dp
441 call scatter(t_Col,t_temp,plan_col3d)
442
443 do i = 1,nlocal
444 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
445 end do
446 endif
447
448 if (do_pot) then
449 ! scatter/gather pot_row into the members of my column
450 call scatter(pot_Row, pot_Temp, plan_row)
451
452 ! scatter/gather pot_local into all other procs
453 ! add resultant to get total pot
454 do i = 1, nlocal
455 pot_local = pot_local + pot_Temp(i)
456 enddo
457
458 pot_Temp = 0.0_DP
459
460 call scatter(pot_Col, pot_Temp, plan_col)
461 do i = 1, nlocal
462 pot_local = pot_local + pot_Temp(i)
463 enddo
464
465 endif
466 #endif
467
468 if (FF_RequiresPostpairCalc() .and. SimRequiresPostpairCalc()) then
469
470 if (FF_uses_RF .and. SimUsesRF()) then
471
472 #ifdef IS_MPI
473 call scatter(rf_Row,rf,plan_row3d)
474 call scatter(rf_Col,rf_Temp,plan_col3d)
475 do i = 1,nlocal
476 rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i)
477 end do
478 #endif
479
480 do i = 1, getNlocal()
481
482 rfpot = 0.0_DP
483 #ifdef IS_MPI
484 me_i = atid_row(i)
485 #else
486 me_i = atid(i)
487 #endif
488 call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
489 if ( is_DP_i ) then
490 call getElementProperty(atypes, me_i, "dipole_moment", mu_i)
491 !! The reaction field needs to include a self contribution
492 !! to the field:
493 call accumulate_self_rf(i, mu_i, u_l)
494 !! Get the reaction field contribution to the
495 !! potential and torques:
496 call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot)
497 #ifdef IS_MPI
498 pot_local = pot_local + rfpot
499 #else
500 pot = pot + rfpot
501
502 #endif
503 endif
504 enddo
505 endif
506 endif
507
508
509 #ifdef IS_MPI
510
511 if (do_pot) then
512 pot = pot + pot_local
513 !! we assume the c code will do the allreduce to get the total potential
514 !! we could do it right here if we needed to...
515 endif
516
517 if (do_stress) then
518 call mpi_allreduce(tau_Temp, tau, 9,mpi_double_precision,mpi_sum, &
519 mpi_comm_world,mpi_err)
520 call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, &
521 mpi_comm_world,mpi_err)
522 endif
523
524 #else
525
526 if (do_stress) then
527 tau = tau_Temp
528 virial = virial_Temp
529 endif
530
531 #endif
532
533 end subroutine do_force_loop
534
535 subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot)
536
537 real( kind = dp ) :: pot
538 real( kind = dp ), dimension(3,getNlocal()) :: u_l
539 real (kind=dp), dimension(9,getNlocal()) :: A
540 real (kind=dp), dimension(3,getNlocal()) :: f
541 real (kind=dp), dimension(3,getNlocal()) :: t
542
543 logical, intent(inout) :: do_pot, do_stress
544 integer, intent(in) :: i, j
545 real ( kind = dp ), intent(inout) :: rijsq
546 real ( kind = dp ) :: r
547 real ( kind = dp ), intent(inout) :: d(3)
548 logical :: is_LJ_i, is_LJ_j
549 logical :: is_DP_i, is_DP_j
550 logical :: is_GB_i, is_GB_j
551 logical :: is_Sticky_i, is_Sticky_j
552 integer :: me_i, me_j
553
554 r = sqrt(rijsq)
555
556 #ifdef IS_MPI
557 if (tagRow(i) .eq. tagColumn(j)) then
558 write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j)
559 endif
560
561 me_i = atid_row(i)
562 me_j = atid_col(j)
563
564 #else
565
566 me_i = atid(i)
567 me_j = atid(j)
568
569 #endif
570
571 if (FF_uses_LJ .and. SimUsesLJ()) then
572 call getElementProperty(atypes, me_i, "is_LJ", is_LJ_i)
573 call getElementProperty(atypes, me_j, "is_LJ", is_LJ_j)
574
575 if ( is_LJ_i .and. is_LJ_j ) &
576 call do_lj_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
577 endif
578
579 if (FF_uses_dipoles .and. SimUsesDipoles()) then
580 call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
581 call getElementProperty(atypes, me_j, "is_DP", is_DP_j)
582
583 if ( is_DP_i .and. is_DP_j ) then
584
585 call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, &
586 do_pot, do_stress)
587 if (FF_uses_RF .and. SimUsesRF()) then
588 call accumulate_rf(i, j, r, u_l)
589 call rf_correct_forces(i, j, d, r, u_l, f, do_stress)
590 endif
591
592 endif
593 endif
594
595 if (FF_uses_Sticky .and. SimUsesSticky()) then
596
597 call getElementProperty(atypes, me_i, "is_Sticky", is_Sticky_i)
598 call getElementProperty(atypes, me_j, "is_Sticky", is_Sticky_j)
599
600 if ( is_Sticky_i .and. is_Sticky_j ) then
601 call do_sticky_pair(i, j, d, r, rijsq, A, pot, f, t, &
602 do_pot, do_stress)
603 endif
604 endif
605
606
607 if (FF_uses_GB .and. SimUsesGB()) then
608
609 call getElementProperty(atypes, me_i, "is_GB", is_GB_i)
610 call getElementProperty(atypes, me_j, "is_GB", is_GB_j)
611
612 if ( is_GB_i .and. is_GB_j ) then
613 call do_gb_pair(i, j, d, r, rijsq, u_l, pot, f, t, &
614 do_pot, do_stress)
615 endif
616 endif
617
618
619
620 end subroutine do_pair
621
622
623 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
624
625 real (kind = dp), dimension(3) :: q_i
626 real (kind = dp), dimension(3) :: q_j
627 real ( kind = dp ), intent(out) :: r_sq
628 real( kind = dp ) :: d(3), scaled(3)
629 integer i
630
631 d(1:3) = q_j(1:3) - q_i(1:3)
632
633 ! Wrap back into periodic box if necessary
634 if ( SimUsesPBC() ) then
635
636 if( .not.boxIsOrthorhombic ) then
637 ! calc the scaled coordinates.
638
639 scaled = matmul(HmatInv, d)
640
641 ! wrap the scaled coordinates
642
643 scaled = scaled - anint(scaled)
644
645
646 ! calc the wrapped real coordinates from the wrapped scaled
647 ! coordinates
648
649 d = matmul(Hmat,scaled)
650
651 else
652 ! calc the scaled coordinates.
653
654 do i = 1, 3
655 scaled(i) = d(i) * HmatInv(i,i)
656
657 ! wrap the scaled coordinates
658
659 scaled(i) = scaled(i) - anint(scaled(i))
660
661 ! calc the wrapped real coordinates from the wrapped scaled
662 ! coordinates
663
664 d(i) = scaled(i)*Hmat(i,i)
665 enddo
666 endif
667
668 endif
669
670 r_sq = dot_product(d,d)
671
672 end subroutine get_interatomic_vector
673
674 subroutine check_initialization(error)
675 integer, intent(out) :: error
676
677 error = 0
678 ! Make sure we are properly initialized.
679 if (.not. do_forces_initialized) then
680 error = -1
681 return
682 endif
683
684 #ifdef IS_MPI
685 if (.not. isMPISimSet()) then
686 write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
687 error = -1
688 return
689 endif
690 #endif
691
692 return
693 end subroutine check_initialization
694
695
696 subroutine zero_work_arrays()
697
698 #ifdef IS_MPI
699
700 q_Row = 0.0_dp
701 q_Col = 0.0_dp
702
703 u_l_Row = 0.0_dp
704 u_l_Col = 0.0_dp
705
706 A_Row = 0.0_dp
707 A_Col = 0.0_dp
708
709 f_Row = 0.0_dp
710 f_Col = 0.0_dp
711 f_Temp = 0.0_dp
712
713 t_Row = 0.0_dp
714 t_Col = 0.0_dp
715 t_Temp = 0.0_dp
716
717 pot_Row = 0.0_dp
718 pot_Col = 0.0_dp
719 pot_Temp = 0.0_dp
720
721 rf_Row = 0.0_dp
722 rf_Col = 0.0_dp
723 rf_Temp = 0.0_dp
724
725 #endif
726
727 rf = 0.0_dp
728 tau_Temp = 0.0_dp
729 virial_Temp = 0.0_dp
730 end subroutine zero_work_arrays
731
732 function skipThisPair(atom1, atom2) result(skip_it)
733 integer, intent(in) :: atom1
734 integer, intent(in), optional :: atom2
735 logical :: skip_it
736 integer :: unique_id_1, unique_id_2
737 integer :: me_i,me_j
738 integer :: i
739
740 skip_it = .false.
741
742 !! there are a number of reasons to skip a pair or a particle
743 !! mostly we do this to exclude atoms who are involved in short
744 !! range interactions (bonds, bends, torsions), but we also need
745 !! to exclude some overcounted interactions that result from
746 !! the parallel decomposition
747
748 #ifdef IS_MPI
749 !! in MPI, we have to look up the unique IDs for each atom
750 unique_id_1 = tagRow(atom1)
751 #else
752 !! in the normal loop, the atom numbers are unique
753 unique_id_1 = atom1
754 #endif
755
756 !! We were called with only one atom, so just check the global exclude
757 !! list for this atom
758 if (.not. present(atom2)) then
759 do i = 1, nExcludes_global
760 if (excludesGlobal(i) == unique_id_1) then
761 skip_it = .true.
762 return
763 end if
764 end do
765 return
766 end if
767
768 #ifdef IS_MPI
769 unique_id_2 = tagColumn(atom2)
770 #else
771 unique_id_2 = atom2
772 #endif
773
774 #ifdef IS_MPI
775 !! this situation should only arise in MPI simulations
776 if (unique_id_1 == unique_id_2) then
777 skip_it = .true.
778 return
779 end if
780
781 !! this prevents us from doing the pair on multiple processors
782 if (unique_id_1 < unique_id_2) then
783 if (mod(unique_id_1 + unique_id_2,2) == 0) then
784 skip_it = .true.
785 return
786 endif
787 else
788 if (mod(unique_id_1 + unique_id_2,2) == 1) then
789 skip_it = .true.
790 return
791 endif
792 endif
793 #endif
794
795 !! the rest of these situations can happen in all simulations:
796 do i = 1, nExcludes_global
797 if ((excludesGlobal(i) == unique_id_1) .or. &
798 (excludesGlobal(i) == unique_id_2)) then
799 skip_it = .true.
800 return
801 endif
802 enddo
803
804 do i = 1, nExcludes_local
805 if (excludesLocal(1,i) == unique_id_1) then
806 if (excludesLocal(2,i) == unique_id_2) then
807 skip_it = .true.
808 return
809 endif
810 else
811 if (excludesLocal(1,i) == unique_id_2) then
812 if (excludesLocal(2,i) == unique_id_1) then
813 skip_it = .true.
814 return
815 endif
816 endif
817 endif
818 end do
819
820 return
821 end function skipThisPair
822
823 function FF_UsesDirectionalAtoms() result(doesit)
824 logical :: doesit
825 doesit = FF_uses_dipoles .or. FF_uses_sticky .or. &
826 FF_uses_GB .or. FF_uses_RF
827 end function FF_UsesDirectionalAtoms
828
829 function FF_RequiresPrepairCalc() result(doesit)
830 logical :: doesit
831 doesit = FF_uses_EAM
832 end function FF_RequiresPrepairCalc
833
834 function FF_RequiresPostpairCalc() result(doesit)
835 logical :: doesit
836 doesit = FF_uses_RF
837 end function FF_RequiresPostpairCalc
838
839 end module do_Forces