ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 486
Committed: Thu Apr 10 16:22:00 2003 UTC (21 years, 3 months ago) by mmeineke
File size: 21261 byte(s)
Log Message:
fixed a bug in symplectic, where presure was only being calculated the first time through.

File Contents

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