ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 626
Committed: Wed Jul 16 21:30:56 2003 UTC (20 years, 11 months ago) by mmeineke
File size: 22195 byte(s)
Log Message:
Changed how cutoffs were handled from C. Now notifyCutoffs in Fortran notifies those who need the information of any changes to cutoffs.

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