ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 1138
Committed: Wed Apr 28 21:39:22 2004 UTC (20 years, 4 months ago) by gezelter
File size: 38526 byte(s)
Log Message:
work on molecular 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 gezelter 1138 !! @version $Id: do_Forces.F90,v 1.52 2004-04-28 21:39:22 gezelter Exp $, $Date: 2004-04-28 21:39:22 $, $Name: not supported by cvs2svn $, $Revision: 1.52 $
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 gezelter 941 use charge_charge
19 mmeineke 377 use reaction_field
20     use gb_pair
21 mmeineke 626 use vector_class
22 chuckv 650 use eam
23 chuckv 657 use status
24 mmeineke 377 #ifdef IS_MPI
25     use mpiSimulation
26     #endif
27    
28     implicit none
29     PRIVATE
30    
31     #define __FORTRAN90
32     #include "fForceField.h"
33    
34 chuckv 900 logical, save :: haveRlist = .false.
35     logical, save :: haveNeighborList = .false.
36 mmeineke 626 logical, save :: havePolicies = .false.
37 chuckv 900 logical, save :: haveSIMvariables = .false.
38     logical, save :: havePropertyMap = .false.
39     logical, save :: haveSaneForceField = .false.
40 mmeineke 377 logical, save :: FF_uses_LJ
41     logical, save :: FF_uses_sticky
42 gezelter 941 logical, save :: FF_uses_charges
43 mmeineke 377 logical, save :: FF_uses_dipoles
44     logical, save :: FF_uses_RF
45     logical, save :: FF_uses_GB
46     logical, save :: FF_uses_EAM
47 chuckv 900 logical, save :: SIM_uses_LJ
48     logical, save :: SIM_uses_sticky
49 gezelter 941 logical, save :: SIM_uses_charges
50 chuckv 900 logical, save :: SIM_uses_dipoles
51     logical, save :: SIM_uses_RF
52     logical, save :: SIM_uses_GB
53     logical, save :: SIM_uses_EAM
54     logical, save :: SIM_requires_postpair_calc
55     logical, save :: SIM_requires_prepair_calc
56     logical, save :: SIM_uses_directional_atoms
57     logical, save :: SIM_uses_PBC
58 gezelter 1138 logical, save :: SIM_uses_molecular_cutoffs
59 mmeineke 377
60 mmeineke 626 real(kind=dp), save :: rlist, rlistsq
61    
62 mmeineke 377 public :: init_FF
63     public :: do_force_loop
64 mmeineke 626 public :: setRlistDF
65 mmeineke 377
66 chuckv 900
67 chuckv 694 #ifdef PROFILE
68 chuckv 883 public :: getforcetime
69     real, save :: forceTime = 0
70     real :: forceTimeInitial, forceTimeFinal
71 mmeineke 891 integer :: nLoops
72 chuckv 694 #endif
73    
74 chuckv 900 type :: Properties
75     logical :: is_lj = .false.
76     logical :: is_sticky = .false.
77     logical :: is_dp = .false.
78     logical :: is_gb = .false.
79     logical :: is_eam = .false.
80 gezelter 941 logical :: is_charge = .false.
81     real(kind=DP) :: charge = 0.0_DP
82 chuckv 900 real(kind=DP) :: dipole_moment = 0.0_DP
83     end type Properties
84 chuckv 895
85 chuckv 900 type(Properties), dimension(:),allocatable :: PropertyMap
86    
87 mmeineke 377 contains
88    
89 mmeineke 626 subroutine setRlistDF( this_rlist )
90    
91     real(kind=dp) :: this_rlist
92    
93     rlist = this_rlist
94     rlistsq = rlist * rlist
95    
96     haveRlist = .true.
97    
98     end subroutine setRlistDF
99    
100 chuckv 900 subroutine createPropertyMap(status)
101     integer :: nAtypes
102     integer :: status
103     integer :: i
104     logical :: thisProperty
105     real (kind=DP) :: thisDPproperty
106    
107     status = 0
108    
109     nAtypes = getSize(atypes)
110    
111     if (nAtypes == 0) then
112     status = -1
113     return
114     end if
115    
116     if (.not. allocated(PropertyMap)) then
117     allocate(PropertyMap(nAtypes))
118     endif
119    
120     do i = 1, nAtypes
121     call getElementProperty(atypes, i, "is_LJ", thisProperty)
122     PropertyMap(i)%is_LJ = thisProperty
123 gezelter 941
124     call getElementProperty(atypes, i, "is_Charge", thisProperty)
125     PropertyMap(i)%is_Charge = thisProperty
126    
127     if (thisProperty) then
128     call getElementProperty(atypes, i, "charge", thisDPproperty)
129     PropertyMap(i)%charge = thisDPproperty
130     endif
131    
132 chuckv 900 call getElementProperty(atypes, i, "is_DP", thisProperty)
133     PropertyMap(i)%is_DP = thisProperty
134    
135     if (thisProperty) then
136     call getElementProperty(atypes, i, "dipole_moment", thisDPproperty)
137     PropertyMap(i)%dipole_moment = thisDPproperty
138     endif
139    
140     call getElementProperty(atypes, i, "is_Sticky", thisProperty)
141     PropertyMap(i)%is_Sticky = thisProperty
142     call getElementProperty(atypes, i, "is_GB", thisProperty)
143     PropertyMap(i)%is_GB = thisProperty
144     call getElementProperty(atypes, i, "is_EAM", thisProperty)
145     PropertyMap(i)%is_EAM = thisProperty
146     end do
147    
148     havePropertyMap = .true.
149    
150     end subroutine createPropertyMap
151    
152     subroutine setSimVariables()
153     SIM_uses_LJ = SimUsesLJ()
154     SIM_uses_sticky = SimUsesSticky()
155 gezelter 941 SIM_uses_charges = SimUsesCharges()
156 chuckv 900 SIM_uses_dipoles = SimUsesDipoles()
157     SIM_uses_RF = SimUsesRF()
158     SIM_uses_GB = SimUsesGB()
159     SIM_uses_EAM = SimUsesEAM()
160     SIM_requires_postpair_calc = SimRequiresPostpairCalc()
161     SIM_requires_prepair_calc = SimRequiresPrepairCalc()
162     SIM_uses_directional_atoms = SimUsesDirectionalAtoms()
163     SIM_uses_PBC = SimUsesPBC()
164 gezelter 1138 SIM_uses_molecular_cutoffs = SimUsesMolecularCutoffs()
165 chuckv 900
166     haveSIMvariables = .true.
167    
168     return
169     end subroutine setSimVariables
170    
171     subroutine doReadyCheck(error)
172     integer, intent(out) :: error
173    
174     integer :: myStatus
175    
176     error = 0
177    
178     if (.not. havePropertyMap) then
179    
180     myStatus = 0
181    
182     call createPropertyMap(myStatus)
183    
184     if (myStatus .ne. 0) then
185     write(default_error, *) 'createPropertyMap failed in do_Forces!'
186     error = -1
187     return
188     endif
189     endif
190    
191     if (.not. haveSIMvariables) then
192     call setSimVariables()
193     endif
194    
195     if (.not. haveRlist) then
196     write(default_error, *) 'rList has not been set in do_Forces!'
197     error = -1
198     return
199     endif
200    
201     if (SIM_uses_LJ .and. FF_uses_LJ) then
202     if (.not. havePolicies) then
203     write(default_error, *) 'LJ mixing Policies have not been set in do_Forces!'
204     error = -1
205     return
206     endif
207     endif
208    
209     if (.not. haveNeighborList) then
210     write(default_error, *) 'neighbor list has not been initialized in do_Forces!'
211     error = -1
212     return
213     end if
214    
215     if (.not. haveSaneForceField) then
216     write(default_error, *) 'Force Field is not sane in do_Forces!'
217     error = -1
218     return
219     end if
220    
221     #ifdef IS_MPI
222     if (.not. isMPISimSet()) then
223     write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
224     error = -1
225     return
226     endif
227     #endif
228     return
229     end subroutine doReadyCheck
230    
231    
232 mmeineke 377 subroutine init_FF(LJMIXPOLICY, use_RF_c, thisStat)
233    
234     integer, intent(in) :: LJMIXPOLICY
235     logical, intent(in) :: use_RF_c
236    
237     integer, intent(out) :: thisStat
238     integer :: my_status, nMatches
239     integer, pointer :: MatchList(:) => null()
240     real(kind=dp) :: rcut, rrf, rt, dielect
241    
242     !! assume things are copacetic, unless they aren't
243     thisStat = 0
244    
245     !! Fortran's version of a cast:
246     FF_uses_RF = use_RF_c
247    
248     !! init_FF is called *after* all of the atom types have been
249     !! defined in atype_module using the new_atype subroutine.
250     !!
251     !! this will scan through the known atypes and figure out what
252     !! interactions are used by the force field.
253    
254     FF_uses_LJ = .false.
255     FF_uses_sticky = .false.
256 gezelter 941 FF_uses_charges = .false.
257 mmeineke 377 FF_uses_dipoles = .false.
258     FF_uses_GB = .false.
259     FF_uses_EAM = .false.
260    
261     call getMatchingElementList(atypes, "is_LJ", .true., nMatches, MatchList)
262     if (nMatches .gt. 0) FF_uses_LJ = .true.
263 gezelter 941
264     call getMatchingElementList(atypes, "is_Charge", .true., nMatches, MatchList)
265 tim 1113 if (nMatches .gt. 0) FF_uses_charges = .true.
266    
267 mmeineke 377 call getMatchingElementList(atypes, "is_DP", .true., nMatches, MatchList)
268     if (nMatches .gt. 0) FF_uses_dipoles = .true.
269    
270     call getMatchingElementList(atypes, "is_Sticky", .true., nMatches, &
271     MatchList)
272     if (nMatches .gt. 0) FF_uses_Sticky = .true.
273    
274     call getMatchingElementList(atypes, "is_GB", .true., nMatches, MatchList)
275     if (nMatches .gt. 0) FF_uses_GB = .true.
276    
277     call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList)
278     if (nMatches .gt. 0) FF_uses_EAM = .true.
279    
280 chuckv 900 !! Assume sanity (for the sake of argument)
281     haveSaneForceField = .true.
282 tim 1132 !!
283     if (FF_uses_charges) then
284     dielect = getDielect()
285     call initialize_charge(dielect)
286     endif
287 chuckv 900
288 tim 1132
289 mmeineke 377 !! check to make sure the FF_uses_RF setting makes sense
290    
291     if (FF_uses_dipoles) then
292     if (FF_uses_RF) then
293     dielect = getDielect()
294 mmeineke 626 call initialize_rf(dielect)
295 mmeineke 377 endif
296     else
297     if (FF_uses_RF) then
298     write(default_error,*) 'Using Reaction Field with no dipoles? Huh?'
299     thisStat = -1
300 chuckv 900 haveSaneForceField = .false.
301 mmeineke 377 return
302     endif
303 mmeineke 626 endif
304 mmeineke 377
305     if (FF_uses_LJ) then
306    
307     select case (LJMIXPOLICY)
308 gezelter 834 case (LB_MIXING_RULE)
309     call init_lj_FF(LB_MIXING_RULE, my_status)
310     case (EXPLICIT_MIXING_RULE)
311     call init_lj_FF(EXPLICIT_MIXING_RULE, my_status)
312 mmeineke 377 case default
313     write(default_error,*) 'unknown LJ Mixing Policy!'
314     thisStat = -1
315 chuckv 900 haveSaneForceField = .false.
316 mmeineke 377 return
317     end select
318     if (my_status /= 0) then
319     thisStat = -1
320 chuckv 900 haveSaneForceField = .false.
321 mmeineke 377 return
322     end if
323 chuckv 900 havePolicies = .true.
324 mmeineke 377 endif
325    
326     if (FF_uses_sticky) then
327     call check_sticky_FF(my_status)
328     if (my_status /= 0) then
329     thisStat = -1
330 chuckv 900 haveSaneForceField = .false.
331 mmeineke 377 return
332     end if
333     endif
334 chuckv 657
335    
336     if (FF_uses_EAM) then
337 chuckv 801 call init_EAM_FF(my_status)
338 chuckv 657 if (my_status /= 0) then
339 chuckv 900 write(default_error, *) "init_EAM_FF returned a bad status"
340 chuckv 657 thisStat = -1
341 chuckv 900 haveSaneForceField = .false.
342 chuckv 657 return
343     end if
344     endif
345    
346 mmeineke 377 if (FF_uses_GB) then
347     call check_gb_pair_FF(my_status)
348     if (my_status .ne. 0) then
349     thisStat = -1
350 chuckv 900 haveSaneForceField = .false.
351 mmeineke 377 return
352     endif
353     endif
354    
355     if (FF_uses_GB .and. FF_uses_LJ) then
356     endif
357 chuckv 900 if (.not. haveNeighborList) then
358 chuckv 480 !! Create neighbor lists
359 chuckv 898 call expandNeighborList(nLocal, my_status)
360 chuckv 480 if (my_Status /= 0) then
361     write(default_error,*) "SimSetup: ExpandNeighborList returned error."
362     thisStat = -1
363     return
364     endif
365 chuckv 900 haveNeighborList = .true.
366 chuckv 480 endif
367 chuckv 657
368 mmeineke 377 end subroutine init_FF
369    
370    
371     !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
372     !------------------------------------------------------------->
373 gezelter 1138 subroutine do_force_loop(q, qcom, mfact, A, u_l, f, t, tau, pot, &
374     do_pot_c, do_stress_c, error)
375 mmeineke 377 !! Position array provided by C, dimensioned by getNlocal
376 chuckv 898 real ( kind = dp ), dimension(3,nLocal) :: q
377 gezelter 1138 !! molecular center-of-mass position array
378     real ( kind = dp ), dimension(3,nLocal) :: qcom
379     !! mass factors used for molecular cutoffs
380     real ( kind = dp ), dimension(3,nLocal) :: mfact
381 mmeineke 377 !! Rotation Matrix for each long range particle in simulation.
382 chuckv 898 real( kind = dp), dimension(9,nLocal) :: A
383 mmeineke 377 !! Unit vectors for dipoles (lab frame)
384 chuckv 898 real( kind = dp ), dimension(3,nLocal) :: u_l
385 mmeineke 377 !! Force array provided by C, dimensioned by getNlocal
386 chuckv 898 real ( kind = dp ), dimension(3,nLocal) :: f
387 mmeineke 377 !! Torsion array provided by C, dimensioned by getNlocal
388 chuckv 898 real( kind = dp ), dimension(3,nLocal) :: t
389 chuckv 895
390 mmeineke 377 !! Stress Tensor
391     real( kind = dp), dimension(9) :: tau
392     real ( kind = dp ) :: pot
393     logical ( kind = 2) :: do_pot_c, do_stress_c
394     logical :: do_pot
395     logical :: do_stress
396 chuckv 439 #ifdef IS_MPI
397 chuckv 441 real( kind = DP ) :: pot_local
398 mmeineke 377 integer :: nrow
399     integer :: ncol
400 chuckv 694 integer :: nprocs
401 mmeineke 377 #endif
402     integer :: natoms
403     logical :: update_nlist
404     integer :: i, j, jbeg, jend, jnab
405     integer :: nlist
406 gezelter 1138 real( kind = DP ) :: rijsq, rcijsq
407     real(kind=dp),dimension(3) :: d, dc
408 mmeineke 377 real(kind=dp) :: rfpot, mu_i, virial
409 chuckv 897 integer :: me_i, me_j
410 mmeineke 377 logical :: is_dp_i
411     integer :: neighborListSize
412     integer :: listerror, error
413     integer :: localError
414 chuckv 897 integer :: propPack_i, propPack_j
415 mmeineke 626
416 gezelter 845 real(kind=dp) :: listSkin = 1.0
417 gezelter 1138
418 mmeineke 377 !! initialize local variables
419 gezelter 1138
420 mmeineke 377 #ifdef IS_MPI
421 chuckv 441 pot_local = 0.0_dp
422 mmeineke 377 nrow = getNrow(plan_row)
423     ncol = getNcol(plan_col)
424     #else
425     natoms = nlocal
426     #endif
427 gezelter 1138
428 chuckv 900 call doReadyCheck(localError)
429 mmeineke 377 if ( localError .ne. 0 ) then
430 chuckv 900 call handleError("do_force_loop", "Not Initialized")
431 mmeineke 377 error = -1
432     return
433     end if
434     call zero_work_arrays()
435 gezelter 1138
436 mmeineke 377 do_pot = do_pot_c
437     do_stress = do_stress_c
438 gezelter 1138
439 mmeineke 377 ! Gather all information needed by all force loops:
440    
441     #ifdef IS_MPI
442 gezelter 1138
443 mmeineke 377 call gather(q,q_Row,plan_row3d)
444     call gather(q,q_Col,plan_col3d)
445 gezelter 1138
446     if (SIM_uses_molecular_cutoffs) then
447     call gather(qcom,qcom_Row,plan_row3d)
448     call gather(qcom,qcom_Col,plan_col3d)
449     endif
450    
451 chuckv 900 if (FF_UsesDirectionalAtoms() .and. SIM_uses_directional_atoms) then
452 mmeineke 377 call gather(u_l,u_l_Row,plan_row3d)
453     call gather(u_l,u_l_Col,plan_col3d)
454    
455     call gather(A,A_Row,plan_row_rotation)
456     call gather(A,A_Col,plan_col_rotation)
457     endif
458    
459     #endif
460 gezelter 1138
461     !! Begin force loop timing:
462 chuckv 694 #ifdef PROFILE
463     call cpu_time(forceTimeInitial)
464     nloops = nloops + 1
465     #endif
466 gezelter 1138
467 chuckv 900 if (FF_RequiresPrepairCalc() .and. SIM_requires_prepair_calc) then
468 mmeineke 377 !! See if we need to update neighbor lists
469 gezelter 1138 if (SIM_uses_molecular_cutoffs) then
470     call checkNeighborList(nlocal, qcom, listSkin, update_nlist)
471     else
472     call checkNeighborList(nlocal, q, listSkin, update_nlist)
473     endif
474 mmeineke 377 !! if_mpi_gather_stuff_for_prepair
475     !! do_prepair_loop_if_needed
476     !! if_mpi_scatter_stuff_from_prepair
477     !! if_mpi_gather_stuff_from_prepair_to_main_loop
478 gezelter 1138
479     !--------------------PREFORCE LOOP----------->>>>>>>>>>>>>>>>>>>>>>>>>>>
480 chuckv 648 #ifdef IS_MPI
481    
482 gezelter 1138 if (update_nlist) then
483 chuckv 648
484 gezelter 1138 !! save current configuration, construct neighbor list,
485     !! and calculate forces
486     if (SIM_uses_molecular_cutoffs) then
487     call saveNeighborList(nlocal, qcom)
488     else
489     call saveNeighborList(nlocal, q)
490     endif
491    
492     neighborListSize = size(list)
493     nlist = 0
494    
495     do i = 1, nrow
496     point(i) = nlist + 1
497 chuckv 648
498 gezelter 1138 prepair_inner: do j = 1, ncol
499 chuckv 648
500 gezelter 1138 if (skipThisPair(i,j)) cycle prepair_inner
501 chuckv 648
502 gezelter 1138 if (SIM_uses_molecular_cutoffs) then
503     call get_interatomic_vector(qcom_Row(:,i), qcom_Col(:,j), &
504     dc, rcijsq)
505     else
506     call get_interatomic_vector(q_Row(:,i), q_Col(:,j), &
507     dc, rcijsq)
508 chuckv 648 endif
509    
510 gezelter 1138 if (rcijsq < rlistsq) then
511    
512     nlist = nlist + 1
513    
514     if (nlist > neighborListSize) then
515     call expandNeighborList(nlocal, listerror)
516     if (listerror /= 0) then
517     error = -1
518     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
519     return
520     end if
521     neighborListSize = size(list)
522     endif
523    
524     list(nlist) = j
525    
526     if (SIM_uses_molecular_cutoffs) then
527     ! since the simulation distances were in molecular COMs:
528     call get_interatomic_vector(q_Row(:,i), q_Col(:,j), &
529     d, rijsq)
530     else
531     d(1:3) = dc(1:3)
532     rijsq = rcijsq
533     endif
534    
535     call do_prepair(i, j, rijsq, d, rcijsq, dc, &
536     do_pot, do_stress, u_l, A, f, t, pot_local)
537     endif
538     enddo prepair_inner
539     enddo
540    
541     point(nrow + 1) = nlist + 1
542    
543     else !! (of update_check)
544    
545     ! use the list to find the neighbors
546     do i = 1, nrow
547     JBEG = POINT(i)
548     JEND = POINT(i+1) - 1
549     ! check thiat molecule i has neighbors
550     if (jbeg .le. jend) then
551    
552     do jnab = jbeg, jend
553     j = list(jnab)
554    
555     call get_interatomic_vector(q_Row(:,i), q_Col(:,j), &
556     d, rijsq)
557     if (SIM_uses_molecular_cutoffs) then
558     call get_interatomic_vector(qcom_Row(:,i),qcom_Col(:,j),&
559     dc, rcijsq)
560     else
561     dc(1:3) = d(1:3)
562     rcijsq = rijsq
563     endif
564    
565     call do_prepair(i, j, rijsq, d, rcijsq, dc, &
566     do_pot, do_stress, u_l, A, f, t, pot_local)
567    
568     enddo
569 chuckv 648 endif
570 gezelter 1138 enddo
571     endif
572 chuckv 648
573     #else
574    
575 gezelter 1138 if (update_nlist) then
576 chuckv 648
577 gezelter 1138 ! save current configuration, contruct neighbor list,
578     ! and calculate forces
579     call saveNeighborList(natoms, q)
580    
581     neighborListSize = size(list)
582    
583     nlist = 0
584    
585     do i = 1, natoms-1
586     point(i) = nlist + 1
587 chuckv 648
588 gezelter 1138 prepair_inner: do j = i+1, natoms
589    
590     if (skipThisPair(i,j)) cycle prepair_inner
591    
592     if (SIM_uses_molecular_cutoffs) then
593     call get_interatomic_vector(qcom(:,i), qcom(:,j), &
594     dc, rcijsq)
595     else
596     call get_interatomic_vector(q(:,i), q(:,j), dc, rcijsq)
597 chuckv 648 endif
598    
599 gezelter 1138 if (rcijsq < rlistsq) then
600    
601     nlist = nlist + 1
602    
603     if (nlist > neighborListSize) then
604     call expandNeighborList(natoms, listerror)
605     if (listerror /= 0) then
606     error = -1
607     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
608     return
609     end if
610     neighborListSize = size(list)
611     endif
612    
613     list(nlist) = j
614    
615    
616     if (SIM_uses_molecular_cutoffs) then
617     ! since the simulation distances were in molecular COMs:
618     call get_interatomic_vector(q(:,i), q(:,j), &
619     d, rijsq)
620     else
621     dc(1:3) = d(1:3)
622     rcijsq = rijsq
623     endif
624    
625     call do_prepair(i, j, rijsq, d, rcijsq, dc, &
626     do_pot, do_stress, u_l, A, f, t, pot)
627    
628     endif
629     enddo prepair_inner
630     enddo
631    
632     point(natoms) = nlist + 1
633    
634     else !! (update)
635    
636     ! use the list to find the neighbors
637     do i = 1, natoms-1
638     JBEG = POINT(i)
639     JEND = POINT(i+1) - 1
640     ! check thiat molecule i has neighbors
641     if (jbeg .le. jend) then
642 chuckv 648
643 gezelter 1138 do jnab = jbeg, jend
644     j = list(jnab)
645    
646     call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
647     if (SIM_uses_molecular_cutoffs) then
648     call get_interatomic_vector(qcom(:,i), qcom(:,j), &
649     dc, rcijsq)
650     else
651     dc(1:3) = d(1:3)
652     rcijsq = rijsq
653     endif
654    
655     call do_prepair(i, j, rijsq, d, rcijsq, dc, &
656     do_pot, do_stress, u_l, A, f, t, pot)
657    
658     enddo
659 chuckv 648 endif
660 gezelter 1138 enddo
661     endif
662 chuckv 648 #endif
663 gezelter 1138 !! Do rest of preforce calculations
664     !! do necessary preforce calculations
665     call do_preforce(nlocal,pot)
666     ! we have already updated the neighbor list set it to false...
667     update_nlist = .false.
668 mmeineke 377 else
669 chuckv 648 !! See if we need to update neighbor lists for non pre-pair
670 mmeineke 626 call checkNeighborList(nlocal, q, listSkin, update_nlist)
671 mmeineke 377 endif
672 gezelter 1138
673     !---------------------------------MAIN Pair LOOP->>>>>>>>>>>>>
674    
675 mmeineke 377 #ifdef IS_MPI
676    
677     if (update_nlist) then
678     !! save current configuration, construct neighbor list,
679     !! and calculate forces
680 mmeineke 459 call saveNeighborList(nlocal, q)
681 mmeineke 377
682     neighborListSize = size(list)
683     nlist = 0
684    
685     do i = 1, nrow
686 gezelter 1138
687 mmeineke 377 point(i) = nlist + 1
688    
689     inner: do j = 1, ncol
690    
691     if (skipThisPair(i,j)) cycle inner
692    
693 gezelter 1138 if (SIM_uses_molecular_cutoffs) then
694     call get_interatomic_vector(qcom_Row(:,i), qcom_Col(:,j), &
695     dc, rcijsq)
696     else
697     call get_interatomic_vector(q_Row(:,i), q_Col(:,j), &
698     dc, rcijsq)
699     endif
700 mmeineke 377
701 gezelter 1138 if (rcijsq < rlistsq) then
702 mmeineke 377
703     nlist = nlist + 1
704    
705     if (nlist > neighborListSize) then
706     call expandNeighborList(nlocal, listerror)
707     if (listerror /= 0) then
708     error = -1
709     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
710     return
711     end if
712     neighborListSize = size(list)
713     endif
714    
715     list(nlist) = j
716 mmeineke 626
717 gezelter 1138 if (SIM_uses_molecular_cutoffs) then
718     call get_interatomic_vector(q_Row(:,i), q_Col(:,j), &
719     d, rijsq)
720     else
721     d(1:3) = dc(1:3)
722     rijsq = rcijsq
723     endif
724    
725     call do_pair(i, j, rijsq, d, rcijsq, dc, mfact, &
726     do_pot, do_stress, u_l, A, f, t, pot_local)
727    
728 mmeineke 377 endif
729     enddo inner
730     enddo
731 gezelter 1138
732 mmeineke 377 point(nrow + 1) = nlist + 1
733    
734     else !! (of update_check)
735 gezelter 1138
736 mmeineke 377 ! use the list to find the neighbors
737     do i = 1, nrow
738     JBEG = POINT(i)
739     JEND = POINT(i+1) - 1
740     ! check thiat molecule i has neighbors
741     if (jbeg .le. jend) then
742    
743     do jnab = jbeg, jend
744     j = list(jnab)
745 gezelter 1138
746 mmeineke 377 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
747 gezelter 1138 if (SIM_uses_molecular_cutoffs) then
748     call get_interatomic_vector(qcom_Row(:,i), qcom_Col(:,j), &
749     dc, rcijsq)
750     else
751     dc(1:3) = d(1:3)
752     rcijsq = rijsq
753     endif
754    
755     call do_pair(i, j, rijsq, d, rcijsq, dc, mfact, &
756     do_pot, do_stress, u_l, A, f, t, pot_local)
757    
758 mmeineke 377 enddo
759     endif
760     enddo
761     endif
762    
763     #else
764    
765     if (update_nlist) then
766 gezelter 1138
767 mmeineke 377 ! save current configuration, contruct neighbor list,
768     ! and calculate forces
769 mmeineke 459 call saveNeighborList(natoms, q)
770 mmeineke 377
771     neighborListSize = size(list)
772 gezelter 1138
773 mmeineke 377 nlist = 0
774    
775     do i = 1, natoms-1
776     point(i) = nlist + 1
777    
778     inner: do j = i+1, natoms
779    
780 chuckv 388 if (skipThisPair(i,j)) cycle inner
781 gezelter 1138
782     if (SIM_uses_molecular_cutoffs) then
783     call get_interatomic_vector(qcom(:,i), qcom(:,j), &
784     dc, rcijsq)
785     else
786     call get_interatomic_vector(q(:,i), q(:,j), &
787     dc, rcijsq)
788     endif
789    
790     if (rcijsq < rlistsq) then
791 mmeineke 377
792     nlist = nlist + 1
793 gezelter 1138
794 mmeineke 377 if (nlist > neighborListSize) then
795     call expandNeighborList(natoms, listerror)
796     if (listerror /= 0) then
797     error = -1
798     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
799     return
800     end if
801     neighborListSize = size(list)
802     endif
803    
804     list(nlist) = j
805    
806 gezelter 1138 if (SIM_uses_molecular_cutoffs) then
807     call get_interatomic_vector(q(:,i), q(:,j), &
808     d, rijsq)
809     else
810     d(1:3) = dc(1:3)
811     rijsq = rcijsq
812     endif
813 mmeineke 626
814 gezelter 1138 call do_pair(i, j, rijsq, d, rcijsq, dc, mfact, &
815     do_pot, do_stress, u_l, A, f, t, pot)
816    
817 mmeineke 377 endif
818     enddo inner
819     enddo
820    
821     point(natoms) = nlist + 1
822    
823     else !! (update)
824    
825     ! use the list to find the neighbors
826     do i = 1, natoms-1
827     JBEG = POINT(i)
828     JEND = POINT(i+1) - 1
829     ! check thiat molecule i has neighbors
830     if (jbeg .le. jend) then
831    
832     do jnab = jbeg, jend
833     j = list(jnab)
834    
835 gezelter 1138
836 mmeineke 377 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
837 gezelter 1138 if (SIM_uses_molecular_cutoffs) then
838     call get_interatomic_vector(qcom(:,i), qcom(:,j), &
839     dc, rcijsq)
840     else
841     dc(1:3) = d(1:3)
842     rcijsq = rijsq
843     endif
844    
845     call do_pair(i, j, rijsq, d, rcijsq, dc, mfact, &
846     do_pot, do_stress, u_l, A, f, t, pot)
847    
848 mmeineke 377 enddo
849     endif
850     enddo
851     endif
852    
853     #endif
854    
855     ! phew, done with main loop.
856 gezelter 1138
857     !! Do timing
858 chuckv 694 #ifdef PROFILE
859     call cpu_time(forceTimeFinal)
860     forceTime = forceTime + forceTimeFinal - forceTimeInitial
861     #endif
862 gezelter 1138
863    
864 mmeineke 377 #ifdef IS_MPI
865     !!distribute forces
866 gezelter 1138
867 chuckv 438 f_temp = 0.0_dp
868     call scatter(f_Row,f_temp,plan_row3d)
869     do i = 1,nlocal
870     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
871     end do
872 gezelter 1138
873 chuckv 438 f_temp = 0.0_dp
874 mmeineke 377 call scatter(f_Col,f_temp,plan_col3d)
875     do i = 1,nlocal
876     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
877     end do
878    
879 chuckv 900 if (FF_UsesDirectionalAtoms() .and. SIM_uses_directional_atoms) then
880 chuckv 438 t_temp = 0.0_dp
881     call scatter(t_Row,t_temp,plan_row3d)
882     do i = 1,nlocal
883     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
884     end do
885     t_temp = 0.0_dp
886 mmeineke 377 call scatter(t_Col,t_temp,plan_col3d)
887    
888     do i = 1,nlocal
889     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
890     end do
891     endif
892    
893     if (do_pot) then
894     ! scatter/gather pot_row into the members of my column
895     call scatter(pot_Row, pot_Temp, plan_row)
896 gezelter 1138
897 mmeineke 377 ! scatter/gather pot_local into all other procs
898     ! add resultant to get total pot
899     do i = 1, nlocal
900     pot_local = pot_local + pot_Temp(i)
901     enddo
902 chuckv 439
903     pot_Temp = 0.0_DP
904 gezelter 1138
905 mmeineke 377 call scatter(pot_Col, pot_Temp, plan_col)
906     do i = 1, nlocal
907     pot_local = pot_local + pot_Temp(i)
908     enddo
909 chuckv 439
910 gezelter 1138 endif
911 mmeineke 377 #endif
912 gezelter 1138
913 chuckv 900 if (FF_RequiresPostpairCalc() .and. SIM_requires_postpair_calc) then
914 mmeineke 377
915 chuckv 900 if (FF_uses_RF .and. SIM_uses_RF) then
916 gezelter 1138
917 mmeineke 377 #ifdef IS_MPI
918     call scatter(rf_Row,rf,plan_row3d)
919     call scatter(rf_Col,rf_Temp,plan_col3d)
920     do i = 1,nlocal
921     rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i)
922     end do
923     #endif
924    
925 chuckv 898 do i = 1, nLocal
926 gezelter 1138
927 mmeineke 377 rfpot = 0.0_DP
928     #ifdef IS_MPI
929     me_i = atid_row(i)
930     #else
931     me_i = atid(i)
932     #endif
933 gezelter 1138
934 chuckv 900 if (PropertyMap(me_i)%is_DP) then
935 gezelter 1138
936 chuckv 900 mu_i = PropertyMap(me_i)%dipole_moment
937 gezelter 1138
938 mmeineke 377 !! The reaction field needs to include a self contribution
939     !! to the field:
940 chuckv 900 call accumulate_self_rf(i, mu_i, u_l)
941 mmeineke 377 !! Get the reaction field contribution to the
942     !! potential and torques:
943     call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot)
944     #ifdef IS_MPI
945     pot_local = pot_local + rfpot
946     #else
947     pot = pot + rfpot
948    
949     #endif
950     endif
951     enddo
952     endif
953     endif
954 gezelter 1138
955    
956 mmeineke 377 #ifdef IS_MPI
957 gezelter 1138
958 mmeineke 377 if (do_pot) then
959 chuckv 441 pot = pot + pot_local
960 mmeineke 377 !! we assume the c code will do the allreduce to get the total potential
961     !! we could do it right here if we needed to...
962     endif
963 gezelter 1138
964 mmeineke 377 if (do_stress) then
965 gezelter 1138 call mpi_allreduce(tau_Temp, tau, 9,mpi_double_precision,mpi_sum, &
966 mmeineke 377 mpi_comm_world,mpi_err)
967 chuckv 470 call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, &
968 mmeineke 377 mpi_comm_world,mpi_err)
969     endif
970 gezelter 1138
971 mmeineke 377 #else
972 gezelter 1138
973 mmeineke 377 if (do_stress) then
974     tau = tau_Temp
975     virial = virial_Temp
976     endif
977 mmeineke 887
978 mmeineke 377 #endif
979 mmeineke 887
980    
981 mmeineke 377 end subroutine do_force_loop
982 gezelter 1138
983     subroutine do_pair(i, j, rijsq, d, rcijsq, dc, mfact, do_pot, do_stress, &
984     u_l, A, f, t, pot)
985 mmeineke 377
986     real( kind = dp ) :: pot
987 gezelter 1138 real( kind = dp ), dimension(nLocal) :: mfact
988 chuckv 898 real( kind = dp ), dimension(3,nLocal) :: u_l
989 gezelter 1138 real( kind = dp ), dimension(9,nLocal) :: A
990     real( kind = dp ), dimension(3,nLocal) :: f
991     real( kind = dp ), dimension(3,nLocal) :: t
992 mmeineke 377
993     logical, intent(inout) :: do_pot, do_stress
994     integer, intent(in) :: i, j
995 gezelter 1138 real ( kind = dp ), intent(inout) :: rijsq, rcijsq
996     real ( kind = dp ) :: r, rc
997     real ( kind = dp ), intent(inout) :: d(3), dc(3)
998 mmeineke 377 integer :: me_i, me_j
999 gezelter 946
1000 mmeineke 377 r = sqrt(rijsq)
1001 gezelter 1138 if (SIM_uses_molecular_cutoffs) then
1002     rc = sqrt(rcijsq)
1003     else
1004     rc = r
1005     endif
1006 mmeineke 377
1007 gezelter 1138
1008 mmeineke 377 #ifdef IS_MPI
1009 gezelter 490 if (tagRow(i) .eq. tagColumn(j)) then
1010     write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j)
1011     endif
1012 mmeineke 377 me_i = atid_row(i)
1013     me_j = atid_col(j)
1014     #else
1015     me_i = atid(i)
1016     me_j = atid(j)
1017     #endif
1018 chuckv 894
1019 chuckv 900 if (FF_uses_LJ .and. SIM_uses_LJ) then
1020 gezelter 946
1021     if ( PropertyMap(me_i)%is_LJ .and. PropertyMap(me_j)%is_LJ ) then
1022     call do_lj_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
1023     endif
1024    
1025 mmeineke 377 endif
1026 gezelter 946
1027     if (FF_uses_charges .and. SIM_uses_charges) then
1028    
1029     if (PropertyMap(me_i)%is_Charge .and. PropertyMap(me_j)%is_Charge) then
1030 gezelter 1138 call do_charge_pair(i, j, d, r, rijsq, dc, rc, rcijsq, mfact, &
1031     pot, f, do_pot, do_stress, SIM_uses_molecular_cutoffs)
1032 gezelter 946 endif
1033    
1034     endif
1035    
1036 chuckv 900 if (FF_uses_dipoles .and. SIM_uses_dipoles) then
1037 mmeineke 377
1038 chuckv 900 if ( PropertyMap(me_i)%is_DP .and. PropertyMap(me_j)%is_DP) then
1039 gezelter 462 call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, &
1040 mmeineke 377 do_pot, do_stress)
1041 chuckv 900 if (FF_uses_RF .and. SIM_uses_RF) then
1042 mmeineke 377 call accumulate_rf(i, j, r, u_l)
1043     call rf_correct_forces(i, j, d, r, u_l, f, do_stress)
1044 gezelter 946 endif
1045 mmeineke 377 endif
1046 gezelter 946
1047 mmeineke 377 endif
1048    
1049 chuckv 900 if (FF_uses_Sticky .and. SIM_uses_sticky) then
1050 mmeineke 377
1051 chuckv 900 if ( PropertyMap(me_i)%is_Sticky .and. PropertyMap(me_j)%is_Sticky) then
1052 mmeineke 377 call do_sticky_pair(i, j, d, r, rijsq, A, pot, f, t, &
1053     do_pot, do_stress)
1054     endif
1055 gezelter 946
1056 mmeineke 377 endif
1057    
1058    
1059 chuckv 900 if (FF_uses_GB .and. SIM_uses_GB) then
1060 mmeineke 377
1061 chuckv 900 if ( PropertyMap(me_i)%is_GB .and. PropertyMap(me_j)%is_GB) then
1062 mmeineke 377 call do_gb_pair(i, j, d, r, rijsq, u_l, pot, f, t, &
1063     do_pot, do_stress)
1064     endif
1065 chuckv 894
1066 mmeineke 377 endif
1067 gezelter 946
1068     if (FF_uses_EAM .and. SIM_uses_EAM) then
1069    
1070     if ( PropertyMap(me_i)%is_EAM .and. PropertyMap(me_j)%is_EAM) then
1071     call do_eam_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
1072     endif
1073    
1074     endif
1075 mmeineke 597
1076 mmeineke 377 end subroutine do_pair
1077    
1078    
1079 chuckv 631
1080 gezelter 1138 subroutine do_prepair(i, j, rijsq, d, rcijsq, dc, &
1081     do_pot, do_stress, u_l, A, f, t, pot)
1082 chuckv 631 real( kind = dp ) :: pot
1083 chuckv 898 real( kind = dp ), dimension(3,nLocal) :: u_l
1084     real (kind=dp), dimension(9,nLocal) :: A
1085     real (kind=dp), dimension(3,nLocal) :: f
1086     real (kind=dp), dimension(3,nLocal) :: t
1087 chuckv 631
1088     logical, intent(inout) :: do_pot, do_stress
1089     integer, intent(in) :: i, j
1090 gezelter 1138 real ( kind = dp ), intent(inout) :: rijsq, rcijsq
1091     real ( kind = dp ) :: r, rc
1092     real ( kind = dp ), intent(inout) :: d(3), dc(3)
1093 chuckv 631
1094     logical :: is_EAM_i, is_EAM_j
1095    
1096     integer :: me_i, me_j
1097    
1098 gezelter 1138
1099     r = sqrt(rijsq)
1100     if (SIM_uses_molecular_cutoffs) then
1101     rc = sqrt(rcijsq)
1102     else
1103     rc = r
1104     endif
1105 chuckv 631
1106 chuckv 669
1107 chuckv 631 #ifdef IS_MPI
1108     if (tagRow(i) .eq. tagColumn(j)) then
1109     write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j)
1110     endif
1111    
1112     me_i = atid_row(i)
1113     me_j = atid_col(j)
1114    
1115     #else
1116    
1117     me_i = atid(i)
1118     me_j = atid(j)
1119    
1120     #endif
1121 chuckv 673
1122 chuckv 900 if (FF_uses_EAM .and. SIM_uses_EAM) then
1123    
1124     if (PropertyMap(me_i)%is_EAM .and. PropertyMap(me_j)%is_EAM) &
1125 chuckv 648 call calc_EAM_prepair_rho(i, j, d, r, rijsq )
1126 chuckv 900
1127 chuckv 631 endif
1128 chuckv 900
1129 chuckv 673 end subroutine do_prepair
1130 chuckv 648
1131    
1132    
1133 chuckv 673
1134 gezelter 1138 subroutine do_preforce(nlocal,pot)
1135     integer :: nlocal
1136     real( kind = dp ) :: pot
1137    
1138     if (FF_uses_EAM .and. SIM_uses_EAM) then
1139     call calc_EAM_preforce_Frho(nlocal,pot)
1140     endif
1141    
1142    
1143     end subroutine do_preforce
1144    
1145    
1146     subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
1147    
1148     real (kind = dp), dimension(3) :: q_i
1149     real (kind = dp), dimension(3) :: q_j
1150     real ( kind = dp ), intent(out) :: r_sq
1151     real( kind = dp ) :: d(3), scaled(3)
1152     integer i
1153    
1154     d(1:3) = q_j(1:3) - q_i(1:3)
1155    
1156     ! Wrap back into periodic box if necessary
1157     if ( SIM_uses_PBC ) then
1158 mmeineke 377
1159 gezelter 1138 if( .not.boxIsOrthorhombic ) then
1160     ! calc the scaled coordinates.
1161    
1162     scaled = matmul(HmatInv, d)
1163    
1164     ! wrap the scaled coordinates
1165    
1166     scaled = scaled - anint(scaled)
1167    
1168    
1169     ! calc the wrapped real coordinates from the wrapped scaled
1170     ! coordinates
1171    
1172     d = matmul(Hmat,scaled)
1173    
1174     else
1175     ! calc the scaled coordinates.
1176    
1177     do i = 1, 3
1178     scaled(i) = d(i) * HmatInv(i,i)
1179    
1180     ! wrap the scaled coordinates
1181    
1182     scaled(i) = scaled(i) - anint(scaled(i))
1183    
1184     ! calc the wrapped real coordinates from the wrapped scaled
1185     ! coordinates
1186    
1187     d(i) = scaled(i)*Hmat(i,i)
1188     enddo
1189     endif
1190    
1191     endif
1192    
1193     r_sq = dot_product(d,d)
1194    
1195     end subroutine get_interatomic_vector
1196 mmeineke 377
1197 gezelter 1138 subroutine zero_work_arrays()
1198    
1199     #ifdef IS_MPI
1200    
1201     q_Row = 0.0_dp
1202     q_Col = 0.0_dp
1203 mmeineke 377
1204 gezelter 1138 qcom_Row = 0.0_dp
1205     qcom_Col = 0.0_dp
1206    
1207     u_l_Row = 0.0_dp
1208     u_l_Col = 0.0_dp
1209    
1210     A_Row = 0.0_dp
1211     A_Col = 0.0_dp
1212    
1213     f_Row = 0.0_dp
1214     f_Col = 0.0_dp
1215     f_Temp = 0.0_dp
1216    
1217     t_Row = 0.0_dp
1218     t_Col = 0.0_dp
1219     t_Temp = 0.0_dp
1220    
1221     pot_Row = 0.0_dp
1222     pot_Col = 0.0_dp
1223     pot_Temp = 0.0_dp
1224    
1225     rf_Row = 0.0_dp
1226     rf_Col = 0.0_dp
1227     rf_Temp = 0.0_dp
1228    
1229 mmeineke 377 #endif
1230 chuckv 673
1231 gezelter 1138 if (FF_uses_EAM .and. SIM_uses_EAM) then
1232     call clean_EAM()
1233     endif
1234    
1235     rf = 0.0_dp
1236     tau_Temp = 0.0_dp
1237     virial_Temp = 0.0_dp
1238     end subroutine zero_work_arrays
1239    
1240     function skipThisPair(atom1, atom2) result(skip_it)
1241     integer, intent(in) :: atom1
1242     integer, intent(in), optional :: atom2
1243     logical :: skip_it
1244     integer :: unique_id_1, unique_id_2
1245     integer :: me_i,me_j
1246     integer :: i
1247    
1248     skip_it = .false.
1249    
1250     !! there are a number of reasons to skip a pair or a particle
1251     !! mostly we do this to exclude atoms who are involved in short
1252     !! range interactions (bonds, bends, torsions), but we also need
1253     !! to exclude some overcounted interactions that result from
1254     !! the parallel decomposition
1255    
1256 mmeineke 377 #ifdef IS_MPI
1257 gezelter 1138 !! in MPI, we have to look up the unique IDs for each atom
1258     unique_id_1 = tagRow(atom1)
1259 mmeineke 377 #else
1260 gezelter 1138 !! in the normal loop, the atom numbers are unique
1261     unique_id_1 = atom1
1262 mmeineke 377 #endif
1263 gezelter 1138
1264     !! We were called with only one atom, so just check the global exclude
1265     !! list for this atom
1266     if (.not. present(atom2)) then
1267     do i = 1, nExcludes_global
1268     if (excludesGlobal(i) == unique_id_1) then
1269     skip_it = .true.
1270     return
1271     end if
1272     end do
1273     return
1274     end if
1275    
1276 mmeineke 377 #ifdef IS_MPI
1277 gezelter 1138 unique_id_2 = tagColumn(atom2)
1278 mmeineke 377 #else
1279 gezelter 1138 unique_id_2 = atom2
1280 mmeineke 377 #endif
1281 gezelter 1138
1282 mmeineke 377 #ifdef IS_MPI
1283 gezelter 1138 !! this situation should only arise in MPI simulations
1284     if (unique_id_1 == unique_id_2) then
1285     skip_it = .true.
1286     return
1287     end if
1288    
1289     !! this prevents us from doing the pair on multiple processors
1290     if (unique_id_1 < unique_id_2) then
1291     if (mod(unique_id_1 + unique_id_2,2) == 0) then
1292     skip_it = .true.
1293     return
1294     endif
1295     else
1296     if (mod(unique_id_1 + unique_id_2,2) == 1) then
1297     skip_it = .true.
1298     return
1299     endif
1300     endif
1301 mmeineke 377 #endif
1302 gezelter 1138
1303     !! the rest of these situations can happen in all simulations:
1304     do i = 1, nExcludes_global
1305     if ((excludesGlobal(i) == unique_id_1) .or. &
1306     (excludesGlobal(i) == unique_id_2)) then
1307     skip_it = .true.
1308     return
1309     endif
1310     enddo
1311    
1312     do i = 1, nExcludes_local
1313     if (excludesLocal(1,i) == unique_id_1) then
1314     if (excludesLocal(2,i) == unique_id_2) then
1315     skip_it = .true.
1316     return
1317     endif
1318     else
1319     if (excludesLocal(1,i) == unique_id_2) then
1320     if (excludesLocal(2,i) == unique_id_1) then
1321     skip_it = .true.
1322     return
1323     endif
1324     endif
1325     endif
1326     end do
1327    
1328     return
1329     end function skipThisPair
1330 chuckv 441
1331 gezelter 1138 function FF_UsesDirectionalAtoms() result(doesit)
1332     logical :: doesit
1333     doesit = FF_uses_dipoles .or. FF_uses_sticky .or. &
1334     FF_uses_GB .or. FF_uses_RF
1335     end function FF_UsesDirectionalAtoms
1336    
1337     function FF_RequiresPrepairCalc() result(doesit)
1338     logical :: doesit
1339     doesit = FF_uses_EAM
1340     end function FF_RequiresPrepairCalc
1341    
1342     function FF_RequiresPostpairCalc() result(doesit)
1343     logical :: doesit
1344     doesit = FF_uses_RF
1345     end function FF_RequiresPostpairCalc
1346    
1347 chuckv 883 #ifdef PROFILE
1348 gezelter 1138 function getforcetime() result(totalforcetime)
1349     real(kind=dp) :: totalforcetime
1350     totalforcetime = forcetime
1351     end function getforcetime
1352 chuckv 883 #endif
1353 gezelter 1138
1354     !! This cleans componets of force arrays belonging only to fortran
1355    
1356 mmeineke 377 end module do_Forces