ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2220
Committed: Thu May 5 14:47:35 2005 UTC (19 years, 2 months ago) by chrisfen
File size: 38399 byte(s)
Log Message:
OOPSE setup for TAP water.  It's not parametrized, but OOPSE will now let me run it...

File Contents

# User Rev Content
1 gezelter 1930 !!
2     !! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3     !!
4     !! The University of Notre Dame grants you ("Licensee") a
5     !! non-exclusive, royalty free, license to use, modify and
6     !! redistribute this software in source and binary code form, provided
7     !! that the following conditions are met:
8     !!
9     !! 1. Acknowledgement of the program authors must be made in any
10     !! publication of scientific results based in part on use of the
11     !! program. An acceptable form of acknowledgement is citation of
12     !! the article in which the program was described (Matthew
13     !! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14     !! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15     !! Parallel Simulation Engine for Molecular Dynamics,"
16     !! J. Comput. Chem. 26, pp. 252-271 (2005))
17     !!
18     !! 2. Redistributions of source code must retain the above copyright
19     !! notice, this list of conditions and the following disclaimer.
20     !!
21     !! 3. Redistributions in binary form must reproduce the above copyright
22     !! notice, this list of conditions and the following disclaimer in the
23     !! documentation and/or other materials provided with the
24     !! distribution.
25     !!
26     !! This software is provided "AS IS," without a warranty of any
27     !! kind. All express or implied conditions, representations and
28     !! warranties, including any implied warranty of merchantability,
29     !! fitness for a particular purpose or non-infringement, are hereby
30     !! excluded. The University of Notre Dame and its licensors shall not
31     !! be liable for any damages suffered by licensee as a result of
32     !! using, modifying or distributing the software or its
33     !! derivatives. In no event will the University of Notre Dame or its
34     !! licensors be liable for any lost revenue, profit or data, or for
35     !! direct, indirect, special, consequential, incidental or punitive
36     !! damages, however caused and regardless of the theory of liability,
37     !! arising out of the use of or inability to use software, even if the
38     !! University of Notre Dame has been advised of the possibility of
39     !! such damages.
40     !!
41    
42 gezelter 1610 !! doForces.F90
43     !! module doForces
44     !! Calculates Long Range forces.
45    
46     !! @author Charles F. Vardeman II
47     !! @author Matthew Meineke
48 chrisfen 2220 !! @version $Id: doForces.F90,v 1.15 2005-05-05 14:47:35 chrisfen Exp $, $Date: 2005-05-05 14:47:35 $, $Name: not supported by cvs2svn $, $Revision: 1.15 $
49 gezelter 1610
50 gezelter 1930
51 gezelter 1610 module doForces
52     use force_globals
53     use simulation
54     use definitions
55     use atype_module
56     use switcheroo
57     use neighborLists
58     use lj
59 gezelter 1930 use sticky
60 gezelter 2085 use electrostatic_module
61 gezelter 1610 use reaction_field
62     use gb_pair
63 chrisfen 1636 use shapes
64 gezelter 1610 use vector_class
65     use eam
66     use status
67     #ifdef IS_MPI
68     use mpiSimulation
69     #endif
70    
71     implicit none
72     PRIVATE
73    
74     #define __FORTRAN90
75     #include "UseTheForce/fSwitchingFunction.h"
76    
77     INTEGER, PARAMETER:: PREPAIR_LOOP = 1
78     INTEGER, PARAMETER:: PAIR_LOOP = 2
79    
80     logical, save :: haveRlist = .false.
81     logical, save :: haveNeighborList = .false.
82     logical, save :: haveSIMvariables = .false.
83     logical, save :: havePropertyMap = .false.
84     logical, save :: haveSaneForceField = .false.
85 gezelter 2204
86 gezelter 1634 logical, save :: FF_uses_DirectionalAtoms
87     logical, save :: FF_uses_LennardJones
88 gezelter 2085 logical, save :: FF_uses_Electrostatics
89     logical, save :: FF_uses_Charges
90     logical, save :: FF_uses_Dipoles
91     logical, save :: FF_uses_Quadrupoles
92 chrisfen 2220 logical, save :: FF_uses_Sticky
93     logical, save :: FF_uses_StickyPower
94 gezelter 1634 logical, save :: FF_uses_GayBerne
95     logical, save :: FF_uses_EAM
96     logical, save :: FF_uses_Shapes
97     logical, save :: FF_uses_FLARB
98 gezelter 1610 logical, save :: FF_uses_RF
99 gezelter 1634
100     logical, save :: SIM_uses_DirectionalAtoms
101     logical, save :: SIM_uses_LennardJones
102     logical, save :: SIM_uses_Electrostatics
103     logical, save :: SIM_uses_Charges
104     logical, save :: SIM_uses_Dipoles
105 gezelter 2085 logical, save :: SIM_uses_Quadrupoles
106 gezelter 1634 logical, save :: SIM_uses_Sticky
107 chrisfen 2220 logical, save :: SIM_uses_StickyPower
108 gezelter 1634 logical, save :: SIM_uses_GayBerne
109     logical, save :: SIM_uses_EAM
110     logical, save :: SIM_uses_Shapes
111     logical, save :: SIM_uses_FLARB
112 gezelter 1610 logical, save :: SIM_uses_RF
113     logical, save :: SIM_requires_postpair_calc
114     logical, save :: SIM_requires_prepair_calc
115     logical, save :: SIM_uses_PBC
116     logical, save :: SIM_uses_molecular_cutoffs
117    
118     real(kind=dp), save :: rlist, rlistsq
119    
120     public :: init_FF
121     public :: do_force_loop
122     public :: setRlistDF
123    
124     #ifdef PROFILE
125     public :: getforcetime
126     real, save :: forceTime = 0
127     real :: forceTimeInitial, forceTimeFinal
128     integer :: nLoops
129     #endif
130    
131     type :: Properties
132 gezelter 1634 logical :: is_Directional = .false.
133     logical :: is_LennardJones = .false.
134     logical :: is_Electrostatic = .false.
135     logical :: is_Charge = .false.
136     logical :: is_Dipole = .false.
137 gezelter 2085 logical :: is_Quadrupole = .false.
138 gezelter 1634 logical :: is_Sticky = .false.
139 chrisfen 2220 logical :: is_StickyPower = .false.
140 gezelter 1634 logical :: is_GayBerne = .false.
141     logical :: is_EAM = .false.
142     logical :: is_Shape = .false.
143     logical :: is_FLARB = .false.
144 gezelter 1610 end type Properties
145    
146     type(Properties), dimension(:),allocatable :: PropertyMap
147    
148     contains
149    
150     subroutine setRlistDF( this_rlist )
151 gezelter 2204
152 gezelter 1610 real(kind=dp) :: this_rlist
153    
154     rlist = this_rlist
155     rlistsq = rlist * rlist
156 gezelter 2204
157 gezelter 1610 haveRlist = .true.
158    
159 gezelter 2204 end subroutine setRlistDF
160 gezelter 1610
161     subroutine createPropertyMap(status)
162     integer :: nAtypes
163     integer :: status
164     integer :: i
165     logical :: thisProperty
166     real (kind=DP) :: thisDPproperty
167    
168     status = 0
169    
170     nAtypes = getSize(atypes)
171    
172     if (nAtypes == 0) then
173     status = -1
174     return
175     end if
176 gezelter 2204
177 gezelter 1610 if (.not. allocated(PropertyMap)) then
178     allocate(PropertyMap(nAtypes))
179     endif
180    
181     do i = 1, nAtypes
182 gezelter 1634 call getElementProperty(atypes, i, "is_Directional", thisProperty)
183     PropertyMap(i)%is_Directional = thisProperty
184 gezelter 1610
185 gezelter 1634 call getElementProperty(atypes, i, "is_LennardJones", thisProperty)
186     PropertyMap(i)%is_LennardJones = thisProperty
187 gezelter 2204
188 gezelter 1634 call getElementProperty(atypes, i, "is_Electrostatic", thisProperty)
189     PropertyMap(i)%is_Electrostatic = thisProperty
190    
191 gezelter 1610 call getElementProperty(atypes, i, "is_Charge", thisProperty)
192     PropertyMap(i)%is_Charge = thisProperty
193 gezelter 2204
194 gezelter 1634 call getElementProperty(atypes, i, "is_Dipole", thisProperty)
195     PropertyMap(i)%is_Dipole = thisProperty
196 gezelter 1610
197 gezelter 2085 call getElementProperty(atypes, i, "is_Quadrupole", thisProperty)
198     PropertyMap(i)%is_Quadrupole = thisProperty
199    
200 gezelter 1634 call getElementProperty(atypes, i, "is_Sticky", thisProperty)
201     PropertyMap(i)%is_Sticky = thisProperty
202 chrisfen 2220
203     call getElementProperty(atypes, i, "is_StickyPower", thisProperty)
204     PropertyMap(i)%is_StickyPower = thisProperty
205 gezelter 1610
206 gezelter 1634 call getElementProperty(atypes, i, "is_GayBerne", thisProperty)
207     PropertyMap(i)%is_GayBerne = thisProperty
208 gezelter 1610
209     call getElementProperty(atypes, i, "is_EAM", thisProperty)
210     PropertyMap(i)%is_EAM = thisProperty
211 gezelter 1634
212     call getElementProperty(atypes, i, "is_Shape", thisProperty)
213     PropertyMap(i)%is_Shape = thisProperty
214    
215     call getElementProperty(atypes, i, "is_FLARB", thisProperty)
216     PropertyMap(i)%is_FLARB = thisProperty
217 gezelter 1610 end do
218    
219     havePropertyMap = .true.
220    
221     end subroutine createPropertyMap
222    
223     subroutine setSimVariables()
224 gezelter 1634 SIM_uses_DirectionalAtoms = SimUsesDirectionalAtoms()
225     SIM_uses_LennardJones = SimUsesLennardJones()
226     SIM_uses_Electrostatics = SimUsesElectrostatics()
227     SIM_uses_Charges = SimUsesCharges()
228     SIM_uses_Dipoles = SimUsesDipoles()
229     SIM_uses_Sticky = SimUsesSticky()
230 chrisfen 2220 SIM_uses_StickyPower = SimUsesStickyPower()
231 gezelter 1634 SIM_uses_GayBerne = SimUsesGayBerne()
232     SIM_uses_EAM = SimUsesEAM()
233     SIM_uses_Shapes = SimUsesShapes()
234     SIM_uses_FLARB = SimUsesFLARB()
235 gezelter 1610 SIM_uses_RF = SimUsesRF()
236     SIM_requires_postpair_calc = SimRequiresPostpairCalc()
237     SIM_requires_prepair_calc = SimRequiresPrepairCalc()
238     SIM_uses_PBC = SimUsesPBC()
239    
240     haveSIMvariables = .true.
241    
242     return
243     end subroutine setSimVariables
244    
245     subroutine doReadyCheck(error)
246     integer, intent(out) :: error
247    
248     integer :: myStatus
249    
250     error = 0
251 gezelter 2204
252 gezelter 1610 if (.not. havePropertyMap) then
253    
254     myStatus = 0
255    
256     call createPropertyMap(myStatus)
257    
258     if (myStatus .ne. 0) then
259     write(default_error, *) 'createPropertyMap failed in doForces!'
260     error = -1
261     return
262     endif
263     endif
264    
265     if (.not. haveSIMvariables) then
266     call setSimVariables()
267     endif
268    
269     if (.not. haveRlist) then
270     write(default_error, *) 'rList has not been set in doForces!'
271     error = -1
272     return
273     endif
274    
275     if (.not. haveNeighborList) then
276     write(default_error, *) 'neighbor list has not been initialized in doForces!'
277     error = -1
278     return
279     end if
280    
281     if (.not. haveSaneForceField) then
282     write(default_error, *) 'Force Field is not sane in doForces!'
283     error = -1
284     return
285     end if
286    
287     #ifdef IS_MPI
288     if (.not. isMPISimSet()) then
289     write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
290     error = -1
291     return
292     endif
293     #endif
294     return
295     end subroutine doReadyCheck
296    
297 gezelter 2204
298 gezelter 1628 subroutine init_FF(use_RF_c, thisStat)
299 gezelter 1610
300     logical, intent(in) :: use_RF_c
301    
302     integer, intent(out) :: thisStat
303     integer :: my_status, nMatches
304     integer, pointer :: MatchList(:) => null()
305     real(kind=dp) :: rcut, rrf, rt, dielect
306    
307     !! assume things are copacetic, unless they aren't
308     thisStat = 0
309    
310     !! Fortran's version of a cast:
311     FF_uses_RF = use_RF_c
312 gezelter 2204
313 gezelter 1610 !! init_FF is called *after* all of the atom types have been
314     !! defined in atype_module using the new_atype subroutine.
315     !!
316     !! this will scan through the known atypes and figure out what
317     !! interactions are used by the force field.
318 gezelter 2204
319 gezelter 1634 FF_uses_DirectionalAtoms = .false.
320     FF_uses_LennardJones = .false.
321 gezelter 2085 FF_uses_Electrostatics = .false.
322 gezelter 1634 FF_uses_Charges = .false.
323     FF_uses_Dipoles = .false.
324     FF_uses_Sticky = .false.
325 chrisfen 2220 FF_uses_StickyPower = .false.
326 gezelter 1634 FF_uses_GayBerne = .false.
327 gezelter 1610 FF_uses_EAM = .false.
328 gezelter 1634 FF_uses_Shapes = .false.
329     FF_uses_FLARB = .false.
330 gezelter 2204
331 gezelter 1634 call getMatchingElementList(atypes, "is_Directional", .true., &
332     nMatches, MatchList)
333     if (nMatches .gt. 0) FF_uses_DirectionalAtoms = .true.
334    
335     call getMatchingElementList(atypes, "is_LennardJones", .true., &
336     nMatches, MatchList)
337     if (nMatches .gt. 0) FF_uses_LennardJones = .true.
338 gezelter 2204
339 gezelter 1634 call getMatchingElementList(atypes, "is_Electrostatic", .true., &
340     nMatches, MatchList)
341     if (nMatches .gt. 0) then
342 gezelter 2085 FF_uses_Electrostatics = .true.
343 gezelter 1634 endif
344    
345     call getMatchingElementList(atypes, "is_Charge", .true., &
346     nMatches, MatchList)
347     if (nMatches .gt. 0) then
348 gezelter 2085 FF_uses_Charges = .true.
349     FF_uses_Electrostatics = .true.
350 gezelter 1634 endif
351 gezelter 2204
352 gezelter 1634 call getMatchingElementList(atypes, "is_Dipole", .true., &
353     nMatches, MatchList)
354     if (nMatches .gt. 0) then
355 gezelter 2085 FF_uses_Dipoles = .true.
356     FF_uses_Electrostatics = .true.
357 gezelter 1634 FF_uses_DirectionalAtoms = .true.
358     endif
359 gezelter 2085
360     call getMatchingElementList(atypes, "is_Quadrupole", .true., &
361     nMatches, MatchList)
362     if (nMatches .gt. 0) then
363     FF_uses_Quadrupoles = .true.
364     FF_uses_Electrostatics = .true.
365     FF_uses_DirectionalAtoms = .true.
366     endif
367 gezelter 2204
368 gezelter 1610 call getMatchingElementList(atypes, "is_Sticky", .true., nMatches, &
369     MatchList)
370 gezelter 1634 if (nMatches .gt. 0) then
371     FF_uses_Sticky = .true.
372     FF_uses_DirectionalAtoms = .true.
373     endif
374 gezelter 2204
375 chrisfen 2220 call getMatchingElementList(atypes, "is_StickyPower", .true., nMatches, &
376     MatchList)
377     if (nMatches .gt. 0) then
378     FF_uses_StickyPower = .true.
379     FF_uses_DirectionalAtoms = .true.
380     endif
381    
382 gezelter 1634 call getMatchingElementList(atypes, "is_GayBerne", .true., &
383     nMatches, MatchList)
384     if (nMatches .gt. 0) then
385     FF_uses_GayBerne = .true.
386     FF_uses_DirectionalAtoms = .true.
387     endif
388 gezelter 2204
389 gezelter 1610 call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList)
390     if (nMatches .gt. 0) FF_uses_EAM = .true.
391 gezelter 2204
392 gezelter 1634 call getMatchingElementList(atypes, "is_Shape", .true., &
393     nMatches, MatchList)
394     if (nMatches .gt. 0) then
395     FF_uses_Shapes = .true.
396     FF_uses_DirectionalAtoms = .true.
397     endif
398    
399     call getMatchingElementList(atypes, "is_FLARB", .true., &
400     nMatches, MatchList)
401     if (nMatches .gt. 0) FF_uses_FLARB = .true.
402    
403 gezelter 1610 !! Assume sanity (for the sake of argument)
404     haveSaneForceField = .true.
405 gezelter 2204
406 gezelter 1610 !! check to make sure the FF_uses_RF setting makes sense
407 gezelter 2204
408 gezelter 1610 if (FF_uses_dipoles) then
409     if (FF_uses_RF) then
410     dielect = getDielect()
411     call initialize_rf(dielect)
412     endif
413     else
414     if (FF_uses_RF) then
415     write(default_error,*) 'Using Reaction Field with no dipoles? Huh?'
416     thisStat = -1
417     haveSaneForceField = .false.
418     return
419     endif
420 gezelter 2204 endif
421 gezelter 1610
422 gezelter 1930 !sticky module does not contain check_sticky_FF anymore
423     !if (FF_uses_sticky) then
424     ! call check_sticky_FF(my_status)
425     ! if (my_status /= 0) then
426     ! thisStat = -1
427     ! haveSaneForceField = .false.
428     ! return
429     ! end if
430     !endif
431 gezelter 1610
432     if (FF_uses_EAM) then
433 gezelter 2204 call init_EAM_FF(my_status)
434 gezelter 1610 if (my_status /= 0) then
435     write(default_error, *) "init_EAM_FF returned a bad status"
436     thisStat = -1
437     haveSaneForceField = .false.
438     return
439     end if
440     endif
441    
442 gezelter 1634 if (FF_uses_GayBerne) then
443 gezelter 1610 call check_gb_pair_FF(my_status)
444     if (my_status .ne. 0) then
445     thisStat = -1
446     haveSaneForceField = .false.
447     return
448     endif
449     endif
450    
451 gezelter 1634 if (FF_uses_GayBerne .and. FF_uses_LennardJones) then
452 gezelter 1610 endif
453 gezelter 2204
454 gezelter 1610 if (.not. haveNeighborList) then
455     !! Create neighbor lists
456     call expandNeighborList(nLocal, my_status)
457     if (my_Status /= 0) then
458     write(default_error,*) "SimSetup: ExpandNeighborList returned error."
459     thisStat = -1
460     return
461     endif
462     haveNeighborList = .true.
463 gezelter 2204 endif
464    
465 gezelter 1610 end subroutine init_FF
466    
467 gezelter 2204
468 gezelter 1610 !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
469     !------------------------------------------------------------->
470 gezelter 1930 subroutine do_force_loop(q, q_group, A, eFrame, f, t, tau, pot, &
471 gezelter 1610 do_pot_c, do_stress_c, error)
472     !! Position array provided by C, dimensioned by getNlocal
473     real ( kind = dp ), dimension(3, nLocal) :: q
474     !! molecular center-of-mass position array
475     real ( kind = dp ), dimension(3, nGroups) :: q_group
476     !! Rotation Matrix for each long range particle in simulation.
477     real( kind = dp), dimension(9, nLocal) :: A
478     !! Unit vectors for dipoles (lab frame)
479 gezelter 1930 real( kind = dp ), dimension(9,nLocal) :: eFrame
480 gezelter 1610 !! Force array provided by C, dimensioned by getNlocal
481     real ( kind = dp ), dimension(3,nLocal) :: f
482     !! Torsion array provided by C, dimensioned by getNlocal
483     real( kind = dp ), dimension(3,nLocal) :: t
484    
485     !! Stress Tensor
486     real( kind = dp), dimension(9) :: tau
487     real ( kind = dp ) :: pot
488     logical ( kind = 2) :: do_pot_c, do_stress_c
489     logical :: do_pot
490     logical :: do_stress
491     logical :: in_switching_region
492     #ifdef IS_MPI
493     real( kind = DP ) :: pot_local
494     integer :: nAtomsInRow
495     integer :: nAtomsInCol
496     integer :: nprocs
497     integer :: nGroupsInRow
498     integer :: nGroupsInCol
499     #endif
500     integer :: natoms
501     logical :: update_nlist
502     integer :: i, j, jstart, jend, jnab
503     integer :: istart, iend
504     integer :: ia, jb, atom1, atom2
505     integer :: nlist
506     real( kind = DP ) :: ratmsq, rgrpsq, rgrp, vpair, vij
507     real( kind = DP ) :: sw, dswdr, swderiv, mf
508     real(kind=dp),dimension(3) :: d_atm, d_grp, fpair, fij
509     real(kind=dp) :: rfpot, mu_i, virial
510     integer :: me_i, me_j, n_in_i, n_in_j
511     logical :: is_dp_i
512     integer :: neighborListSize
513     integer :: listerror, error
514     integer :: localError
515     integer :: propPack_i, propPack_j
516     integer :: loopStart, loopEnd, loop
517    
518     real(kind=dp) :: listSkin = 1.0
519 gezelter 2204
520 gezelter 1610 !! initialize local variables
521 gezelter 2204
522 gezelter 1610 #ifdef IS_MPI
523     pot_local = 0.0_dp
524     nAtomsInRow = getNatomsInRow(plan_atom_row)
525     nAtomsInCol = getNatomsInCol(plan_atom_col)
526     nGroupsInRow = getNgroupsInRow(plan_group_row)
527     nGroupsInCol = getNgroupsInCol(plan_group_col)
528     #else
529     natoms = nlocal
530     #endif
531 gezelter 2204
532 gezelter 1610 call doReadyCheck(localError)
533     if ( localError .ne. 0 ) then
534     call handleError("do_force_loop", "Not Initialized")
535     error = -1
536     return
537     end if
538     call zero_work_arrays()
539 gezelter 2204
540 gezelter 1610 do_pot = do_pot_c
541     do_stress = do_stress_c
542 gezelter 2204
543 gezelter 1610 ! Gather all information needed by all force loops:
544 gezelter 2204
545 gezelter 1610 #ifdef IS_MPI
546 gezelter 2204
547 gezelter 1610 call gather(q, q_Row, plan_atom_row_3d)
548     call gather(q, q_Col, plan_atom_col_3d)
549    
550     call gather(q_group, q_group_Row, plan_group_row_3d)
551     call gather(q_group, q_group_Col, plan_group_col_3d)
552 gezelter 2204
553 gezelter 1634 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
554 gezelter 1930 call gather(eFrame, eFrame_Row, plan_atom_row_rotation)
555     call gather(eFrame, eFrame_Col, plan_atom_col_rotation)
556 gezelter 2204
557 gezelter 1610 call gather(A, A_Row, plan_atom_row_rotation)
558     call gather(A, A_Col, plan_atom_col_rotation)
559     endif
560 gezelter 2204
561 gezelter 1610 #endif
562 gezelter 2204
563 gezelter 1610 !! Begin force loop timing:
564     #ifdef PROFILE
565     call cpu_time(forceTimeInitial)
566     nloops = nloops + 1
567     #endif
568 gezelter 2204
569 gezelter 1610 loopEnd = PAIR_LOOP
570     if (FF_RequiresPrepairCalc() .and. SIM_requires_prepair_calc) then
571     loopStart = PREPAIR_LOOP
572     else
573     loopStart = PAIR_LOOP
574     endif
575    
576     do loop = loopStart, loopEnd
577    
578     ! See if we need to update neighbor lists
579     ! (but only on the first time through):
580     if (loop .eq. loopStart) then
581     #ifdef IS_MPI
582     call checkNeighborList(nGroupsInRow, q_group_row, listSkin, &
583 gezelter 2204 update_nlist)
584 gezelter 1610 #else
585     call checkNeighborList(nGroups, q_group, listSkin, &
586 gezelter 2204 update_nlist)
587 gezelter 1610 #endif
588     endif
589 gezelter 2204
590 gezelter 1610 if (update_nlist) then
591     !! save current configuration and construct neighbor list
592     #ifdef IS_MPI
593     call saveNeighborList(nGroupsInRow, q_group_row)
594     #else
595     call saveNeighborList(nGroups, q_group)
596     #endif
597     neighborListSize = size(list)
598     nlist = 0
599     endif
600 gezelter 2204
601 gezelter 1610 istart = 1
602     #ifdef IS_MPI
603     iend = nGroupsInRow
604     #else
605     iend = nGroups - 1
606     #endif
607     outer: do i = istart, iend
608    
609     if (update_nlist) point(i) = nlist + 1
610 gezelter 2204
611 gezelter 1610 n_in_i = groupStartRow(i+1) - groupStartRow(i)
612 gezelter 2204
613 gezelter 1610 if (update_nlist) then
614     #ifdef IS_MPI
615     jstart = 1
616     jend = nGroupsInCol
617     #else
618     jstart = i+1
619     jend = nGroups
620     #endif
621     else
622     jstart = point(i)
623     jend = point(i+1) - 1
624     ! make sure group i has neighbors
625     if (jstart .gt. jend) cycle outer
626     endif
627 gezelter 2204
628 gezelter 1610 do jnab = jstart, jend
629     if (update_nlist) then
630     j = jnab
631     else
632     j = list(jnab)
633     endif
634    
635     #ifdef IS_MPI
636     call get_interatomic_vector(q_group_Row(:,i), &
637     q_group_Col(:,j), d_grp, rgrpsq)
638     #else
639     call get_interatomic_vector(q_group(:,i), &
640     q_group(:,j), d_grp, rgrpsq)
641     #endif
642    
643     if (rgrpsq < rlistsq) then
644     if (update_nlist) then
645     nlist = nlist + 1
646 gezelter 2204
647 gezelter 1610 if (nlist > neighborListSize) then
648     #ifdef IS_MPI
649     call expandNeighborList(nGroupsInRow, listerror)
650     #else
651     call expandNeighborList(nGroups, listerror)
652     #endif
653     if (listerror /= 0) then
654     error = -1
655     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
656     return
657     end if
658     neighborListSize = size(list)
659     endif
660 gezelter 2204
661 gezelter 1610 list(nlist) = j
662     endif
663 gezelter 2204
664 gezelter 1610 if (loop .eq. PAIR_LOOP) then
665     vij = 0.0d0
666     fij(1:3) = 0.0d0
667     endif
668 gezelter 2204
669 gezelter 1610 call get_switch(rgrpsq, sw, dswdr, rgrp, group_switch, &
670     in_switching_region)
671 gezelter 2204
672 gezelter 1610 n_in_j = groupStartCol(j+1) - groupStartCol(j)
673 gezelter 2204
674 gezelter 1610 do ia = groupStartRow(i), groupStartRow(i+1)-1
675 gezelter 2204
676 gezelter 1610 atom1 = groupListRow(ia)
677 gezelter 2204
678 gezelter 1610 inner: do jb = groupStartCol(j), groupStartCol(j+1)-1
679 gezelter 2204
680 gezelter 1610 atom2 = groupListCol(jb)
681 gezelter 2204
682 gezelter 1610 if (skipThisPair(atom1, atom2)) cycle inner
683    
684     if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
685     d_atm(1:3) = d_grp(1:3)
686     ratmsq = rgrpsq
687     else
688     #ifdef IS_MPI
689     call get_interatomic_vector(q_Row(:,atom1), &
690     q_Col(:,atom2), d_atm, ratmsq)
691     #else
692     call get_interatomic_vector(q(:,atom1), &
693     q(:,atom2), d_atm, ratmsq)
694     #endif
695     endif
696    
697     if (loop .eq. PREPAIR_LOOP) then
698     #ifdef IS_MPI
699     call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
700     rgrpsq, d_grp, do_pot, do_stress, &
701 gezelter 1930 eFrame, A, f, t, pot_local)
702 gezelter 1610 #else
703     call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
704     rgrpsq, d_grp, do_pot, do_stress, &
705 gezelter 1930 eFrame, A, f, t, pot)
706 gezelter 1610 #endif
707     else
708     #ifdef IS_MPI
709     call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
710     do_pot, &
711 gezelter 1930 eFrame, A, f, t, pot_local, vpair, fpair)
712 gezelter 1610 #else
713     call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
714     do_pot, &
715 gezelter 1930 eFrame, A, f, t, pot, vpair, fpair)
716 gezelter 1610 #endif
717    
718     vij = vij + vpair
719     fij(1:3) = fij(1:3) + fpair(1:3)
720     endif
721     enddo inner
722     enddo
723 gezelter 2204
724 gezelter 1610 if (loop .eq. PAIR_LOOP) then
725     if (in_switching_region) then
726     swderiv = vij*dswdr/rgrp
727     fij(1) = fij(1) + swderiv*d_grp(1)
728     fij(2) = fij(2) + swderiv*d_grp(2)
729     fij(3) = fij(3) + swderiv*d_grp(3)
730 gezelter 2204
731 gezelter 1610 do ia=groupStartRow(i), groupStartRow(i+1)-1
732     atom1=groupListRow(ia)
733     mf = mfactRow(atom1)
734     #ifdef IS_MPI
735     f_Row(1,atom1) = f_Row(1,atom1) + swderiv*d_grp(1)*mf
736     f_Row(2,atom1) = f_Row(2,atom1) + swderiv*d_grp(2)*mf
737     f_Row(3,atom1) = f_Row(3,atom1) + swderiv*d_grp(3)*mf
738     #else
739     f(1,atom1) = f(1,atom1) + swderiv*d_grp(1)*mf
740     f(2,atom1) = f(2,atom1) + swderiv*d_grp(2)*mf
741     f(3,atom1) = f(3,atom1) + swderiv*d_grp(3)*mf
742     #endif
743     enddo
744 gezelter 2204
745 gezelter 1610 do jb=groupStartCol(j), groupStartCol(j+1)-1
746     atom2=groupListCol(jb)
747     mf = mfactCol(atom2)
748     #ifdef IS_MPI
749     f_Col(1,atom2) = f_Col(1,atom2) - swderiv*d_grp(1)*mf
750     f_Col(2,atom2) = f_Col(2,atom2) - swderiv*d_grp(2)*mf
751     f_Col(3,atom2) = f_Col(3,atom2) - swderiv*d_grp(3)*mf
752     #else
753     f(1,atom2) = f(1,atom2) - swderiv*d_grp(1)*mf
754     f(2,atom2) = f(2,atom2) - swderiv*d_grp(2)*mf
755     f(3,atom2) = f(3,atom2) - swderiv*d_grp(3)*mf
756     #endif
757     enddo
758     endif
759 gezelter 2204
760 gezelter 1610 if (do_stress) call add_stress_tensor(d_grp, fij)
761     endif
762     end if
763     enddo
764     enddo outer
765 gezelter 2204
766 gezelter 1610 if (update_nlist) then
767     #ifdef IS_MPI
768     point(nGroupsInRow + 1) = nlist + 1
769     #else
770     point(nGroups) = nlist + 1
771     #endif
772     if (loop .eq. PREPAIR_LOOP) then
773     ! we just did the neighbor list update on the first
774     ! pass, so we don't need to do it
775     ! again on the second pass
776     update_nlist = .false.
777     endif
778     endif
779 gezelter 2204
780 gezelter 1610 if (loop .eq. PREPAIR_LOOP) then
781     call do_preforce(nlocal, pot)
782     endif
783 gezelter 2204
784 gezelter 1610 enddo
785 gezelter 2204
786 gezelter 1610 !! Do timing
787     #ifdef PROFILE
788     call cpu_time(forceTimeFinal)
789     forceTime = forceTime + forceTimeFinal - forceTimeInitial
790     #endif
791 gezelter 2204
792 gezelter 1610 #ifdef IS_MPI
793     !!distribute forces
794 gezelter 2204
795 gezelter 1610 f_temp = 0.0_dp
796     call scatter(f_Row,f_temp,plan_atom_row_3d)
797     do i = 1,nlocal
798     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
799     end do
800 gezelter 2204
801 gezelter 1610 f_temp = 0.0_dp
802     call scatter(f_Col,f_temp,plan_atom_col_3d)
803     do i = 1,nlocal
804     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
805     end do
806 gezelter 2204
807 gezelter 1634 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
808 gezelter 1610 t_temp = 0.0_dp
809     call scatter(t_Row,t_temp,plan_atom_row_3d)
810     do i = 1,nlocal
811     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
812     end do
813     t_temp = 0.0_dp
814     call scatter(t_Col,t_temp,plan_atom_col_3d)
815 gezelter 2204
816 gezelter 1610 do i = 1,nlocal
817     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
818     end do
819     endif
820 gezelter 2204
821 gezelter 1610 if (do_pot) then
822     ! scatter/gather pot_row into the members of my column
823     call scatter(pot_Row, pot_Temp, plan_atom_row)
824 gezelter 2204
825 gezelter 1610 ! scatter/gather pot_local into all other procs
826     ! add resultant to get total pot
827     do i = 1, nlocal
828     pot_local = pot_local + pot_Temp(i)
829     enddo
830 gezelter 2204
831 gezelter 1610 pot_Temp = 0.0_DP
832 gezelter 2204
833 gezelter 1610 call scatter(pot_Col, pot_Temp, plan_atom_col)
834     do i = 1, nlocal
835     pot_local = pot_local + pot_Temp(i)
836     enddo
837 gezelter 2204
838 gezelter 1610 endif
839     #endif
840 gezelter 2204
841 gezelter 1610 if (FF_RequiresPostpairCalc() .and. SIM_requires_postpair_calc) then
842 gezelter 2204
843 gezelter 1610 if (FF_uses_RF .and. SIM_uses_RF) then
844 gezelter 2204
845 gezelter 1610 #ifdef IS_MPI
846     call scatter(rf_Row,rf,plan_atom_row_3d)
847     call scatter(rf_Col,rf_Temp,plan_atom_col_3d)
848     do i = 1,nlocal
849     rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i)
850     end do
851     #endif
852 gezelter 2204
853 gezelter 1610 do i = 1, nLocal
854 gezelter 2204
855 gezelter 1610 rfpot = 0.0_DP
856     #ifdef IS_MPI
857     me_i = atid_row(i)
858     #else
859     me_i = atid(i)
860     #endif
861 gezelter 2204
862 gezelter 1634 if (PropertyMap(me_i)%is_Dipole) then
863 gezelter 2204
864 gezelter 1634 mu_i = getDipoleMoment(me_i)
865 gezelter 2204
866 gezelter 1610 !! The reaction field needs to include a self contribution
867     !! to the field:
868 gezelter 1930 call accumulate_self_rf(i, mu_i, eFrame)
869 gezelter 1610 !! Get the reaction field contribution to the
870     !! potential and torques:
871 gezelter 1930 call reaction_field_final(i, mu_i, eFrame, rfpot, t, do_pot)
872 gezelter 1610 #ifdef IS_MPI
873     pot_local = pot_local + rfpot
874     #else
875     pot = pot + rfpot
876 gezelter 2204
877 gezelter 1610 #endif
878 gezelter 2204 endif
879 gezelter 1610 enddo
880     endif
881     endif
882 gezelter 2204
883    
884 gezelter 1610 #ifdef IS_MPI
885 gezelter 2204
886 gezelter 1610 if (do_pot) then
887     pot = pot + pot_local
888     !! we assume the c code will do the allreduce to get the total potential
889     !! we could do it right here if we needed to...
890     endif
891 gezelter 2204
892 gezelter 1610 if (do_stress) then
893     call mpi_allreduce(tau_Temp, tau, 9,mpi_double_precision,mpi_sum, &
894     mpi_comm_world,mpi_err)
895     call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, &
896     mpi_comm_world,mpi_err)
897     endif
898 gezelter 2204
899 gezelter 1610 #else
900 gezelter 2204
901 gezelter 1610 if (do_stress) then
902     tau = tau_Temp
903     virial = virial_Temp
904     endif
905 gezelter 2204
906 gezelter 1610 #endif
907 gezelter 2204
908 gezelter 1610 end subroutine do_force_loop
909 gezelter 2204
910 gezelter 1610 subroutine do_pair(i, j, rijsq, d, sw, do_pot, &
911 gezelter 1930 eFrame, A, f, t, pot, vpair, fpair)
912 gezelter 1610
913     real( kind = dp ) :: pot, vpair, sw
914     real( kind = dp ), dimension(3) :: fpair
915     real( kind = dp ), dimension(nLocal) :: mfact
916 gezelter 1930 real( kind = dp ), dimension(9,nLocal) :: eFrame
917 gezelter 1610 real( kind = dp ), dimension(9,nLocal) :: A
918     real( kind = dp ), dimension(3,nLocal) :: f
919     real( kind = dp ), dimension(3,nLocal) :: t
920    
921     logical, intent(inout) :: do_pot
922     integer, intent(in) :: i, j
923     real ( kind = dp ), intent(inout) :: rijsq
924     real ( kind = dp ) :: r
925     real ( kind = dp ), intent(inout) :: d(3)
926     integer :: me_i, me_j
927    
928     r = sqrt(rijsq)
929     vpair = 0.0d0
930     fpair(1:3) = 0.0d0
931    
932     #ifdef IS_MPI
933     me_i = atid_row(i)
934     me_j = atid_col(j)
935     #else
936     me_i = atid(i)
937     me_j = atid(j)
938     #endif
939 gezelter 1706
940 gezelter 2204 ! write(*,*) i, j, me_i, me_j
941    
942 gezelter 1634 if (FF_uses_LennardJones .and. SIM_uses_LennardJones) then
943 gezelter 2204
944 gezelter 1634 if ( PropertyMap(me_i)%is_LennardJones .and. &
945     PropertyMap(me_j)%is_LennardJones ) then
946 gezelter 1610 call do_lj_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, f, do_pot)
947     endif
948 gezelter 2204
949 gezelter 1610 endif
950 gezelter 2204
951 gezelter 2085 if (FF_uses_Electrostatics .and. SIM_uses_Electrostatics) then
952 gezelter 2204
953 gezelter 2085 if (PropertyMap(me_i)%is_Electrostatic .and. &
954     PropertyMap(me_j)%is_Electrostatic) then
955     call doElectrostaticPair(i, j, d, r, rijsq, sw, vpair, fpair, &
956     pot, eFrame, f, t, do_pot)
957 gezelter 1610 endif
958 gezelter 2204
959 gezelter 2085 if (FF_uses_dipoles .and. SIM_uses_dipoles) then
960     if ( PropertyMap(me_i)%is_Dipole .and. &
961     PropertyMap(me_j)%is_Dipole) then
962     if (FF_uses_RF .and. SIM_uses_RF) then
963     call accumulate_rf(i, j, r, eFrame, sw)
964     call rf_correct_forces(i, j, d, r, eFrame, sw, f, fpair)
965     endif
966 gezelter 1634 endif
967 gezelter 1610 endif
968     endif
969    
970 gezelter 2085
971 gezelter 1610 if (FF_uses_Sticky .and. SIM_uses_sticky) then
972    
973     if ( PropertyMap(me_i)%is_Sticky .and. PropertyMap(me_j)%is_Sticky) then
974 gezelter 1634 call do_sticky_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
975     pot, A, f, t, do_pot)
976 gezelter 1610 endif
977 gezelter 2204
978 gezelter 1610 endif
979    
980 chrisfen 2220 if (FF_uses_StickyPower .and. SIM_uses_stickypower) then
981     if ( PropertyMap(me_i)%is_StickyPower .and. &
982     PropertyMap(me_j)%is_StickyPower) then
983     call do_sticky_power_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
984     pot, A, f, t, do_pot)
985     endif
986     endif
987    
988 gezelter 1634 if (FF_uses_GayBerne .and. SIM_uses_GayBerne) then
989 gezelter 2204
990 gezelter 1634 if ( PropertyMap(me_i)%is_GayBerne .and. &
991     PropertyMap(me_j)%is_GayBerne) then
992     call do_gb_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
993 gezelter 1930 pot, A, f, t, do_pot)
994 gezelter 1610 endif
995 gezelter 2204
996 gezelter 1610 endif
997 gezelter 2204
998 gezelter 1610 if (FF_uses_EAM .and. SIM_uses_EAM) then
999 gezelter 2204
1000 gezelter 1610 if ( PropertyMap(me_i)%is_EAM .and. PropertyMap(me_j)%is_EAM) then
1001     call do_eam_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, f, &
1002     do_pot)
1003     endif
1004 gezelter 2204
1005 gezelter 1610 endif
1006 gezelter 1634
1007 gezelter 1706
1008 gezelter 2204 ! write(*,*) PropertyMap(me_i)%is_Shape,PropertyMap(me_j)%is_Shape
1009 gezelter 1706
1010 gezelter 1634 if (FF_uses_Shapes .and. SIM_uses_Shapes) then
1011     if ( PropertyMap(me_i)%is_Shape .and. &
1012     PropertyMap(me_j)%is_Shape ) then
1013     call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1014 gezelter 1650 pot, A, f, t, do_pot)
1015 gezelter 1634 endif
1016 chrisfen 2211 if ( (PropertyMap(me_i)%is_Shape .and. &
1017     PropertyMap(me_j)%is_LennardJones) .or. &
1018     (PropertyMap(me_i)%is_LennardJones .and. &
1019     PropertyMap(me_j)%is_Shape) ) then
1020     call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1021     pot, A, f, t, do_pot)
1022     endif
1023 gezelter 1634 endif
1024 gezelter 2204
1025 gezelter 1610 end subroutine do_pair
1026    
1027     subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, &
1028 gezelter 1930 do_pot, do_stress, eFrame, A, f, t, pot)
1029 gezelter 1610
1030 gezelter 2204 real( kind = dp ) :: pot, sw
1031     real( kind = dp ), dimension(9,nLocal) :: eFrame
1032     real (kind=dp), dimension(9,nLocal) :: A
1033     real (kind=dp), dimension(3,nLocal) :: f
1034     real (kind=dp), dimension(3,nLocal) :: t
1035 gezelter 1610
1036 gezelter 2204 logical, intent(inout) :: do_pot, do_stress
1037     integer, intent(in) :: i, j
1038     real ( kind = dp ), intent(inout) :: rijsq, rcijsq
1039     real ( kind = dp ) :: r, rc
1040     real ( kind = dp ), intent(inout) :: d(3), dc(3)
1041    
1042     logical :: is_EAM_i, is_EAM_j
1043    
1044     integer :: me_i, me_j
1045    
1046    
1047 gezelter 1610 r = sqrt(rijsq)
1048     if (SIM_uses_molecular_cutoffs) then
1049     rc = sqrt(rcijsq)
1050     else
1051     rc = r
1052     endif
1053    
1054 gezelter 2204
1055 gezelter 1610 #ifdef IS_MPI
1056 gezelter 2204 me_i = atid_row(i)
1057     me_j = atid_col(j)
1058 gezelter 1610 #else
1059 gezelter 2204 me_i = atid(i)
1060     me_j = atid(j)
1061 gezelter 1610 #endif
1062 gezelter 2204
1063     if (FF_uses_EAM .and. SIM_uses_EAM) then
1064    
1065     if (PropertyMap(me_i)%is_EAM .and. PropertyMap(me_j)%is_EAM) &
1066     call calc_EAM_prepair_rho(i, j, d, r, rijsq )
1067    
1068     endif
1069    
1070     end subroutine do_prepair
1071    
1072    
1073     subroutine do_preforce(nlocal,pot)
1074     integer :: nlocal
1075     real( kind = dp ) :: pot
1076    
1077     if (FF_uses_EAM .and. SIM_uses_EAM) then
1078     call calc_EAM_preforce_Frho(nlocal,pot)
1079     endif
1080    
1081    
1082     end subroutine do_preforce
1083    
1084    
1085     subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
1086    
1087     real (kind = dp), dimension(3) :: q_i
1088     real (kind = dp), dimension(3) :: q_j
1089     real ( kind = dp ), intent(out) :: r_sq
1090     real( kind = dp ) :: d(3), scaled(3)
1091     integer i
1092    
1093     d(1:3) = q_j(1:3) - q_i(1:3)
1094    
1095     ! Wrap back into periodic box if necessary
1096     if ( SIM_uses_PBC ) then
1097    
1098     if( .not.boxIsOrthorhombic ) then
1099     ! calc the scaled coordinates.
1100    
1101     scaled = matmul(HmatInv, d)
1102    
1103     ! wrap the scaled coordinates
1104    
1105     scaled = scaled - anint(scaled)
1106    
1107    
1108     ! calc the wrapped real coordinates from the wrapped scaled
1109     ! coordinates
1110    
1111     d = matmul(Hmat,scaled)
1112    
1113     else
1114     ! calc the scaled coordinates.
1115    
1116     do i = 1, 3
1117     scaled(i) = d(i) * HmatInv(i,i)
1118    
1119     ! wrap the scaled coordinates
1120    
1121     scaled(i) = scaled(i) - anint(scaled(i))
1122    
1123     ! calc the wrapped real coordinates from the wrapped scaled
1124     ! coordinates
1125    
1126     d(i) = scaled(i)*Hmat(i,i)
1127     enddo
1128     endif
1129    
1130     endif
1131    
1132     r_sq = dot_product(d,d)
1133    
1134     end subroutine get_interatomic_vector
1135    
1136     subroutine zero_work_arrays()
1137    
1138 gezelter 1610 #ifdef IS_MPI
1139    
1140 gezelter 2204 q_Row = 0.0_dp
1141     q_Col = 0.0_dp
1142    
1143     q_group_Row = 0.0_dp
1144     q_group_Col = 0.0_dp
1145    
1146     eFrame_Row = 0.0_dp
1147     eFrame_Col = 0.0_dp
1148    
1149     A_Row = 0.0_dp
1150     A_Col = 0.0_dp
1151    
1152     f_Row = 0.0_dp
1153     f_Col = 0.0_dp
1154     f_Temp = 0.0_dp
1155    
1156     t_Row = 0.0_dp
1157     t_Col = 0.0_dp
1158     t_Temp = 0.0_dp
1159    
1160     pot_Row = 0.0_dp
1161     pot_Col = 0.0_dp
1162     pot_Temp = 0.0_dp
1163    
1164     rf_Row = 0.0_dp
1165     rf_Col = 0.0_dp
1166     rf_Temp = 0.0_dp
1167    
1168 gezelter 1610 #endif
1169 gezelter 2204
1170     if (FF_uses_EAM .and. SIM_uses_EAM) then
1171     call clean_EAM()
1172     endif
1173    
1174     rf = 0.0_dp
1175     tau_Temp = 0.0_dp
1176     virial_Temp = 0.0_dp
1177     end subroutine zero_work_arrays
1178    
1179     function skipThisPair(atom1, atom2) result(skip_it)
1180     integer, intent(in) :: atom1
1181     integer, intent(in), optional :: atom2
1182     logical :: skip_it
1183     integer :: unique_id_1, unique_id_2
1184     integer :: me_i,me_j
1185     integer :: i
1186    
1187     skip_it = .false.
1188    
1189     !! there are a number of reasons to skip a pair or a particle
1190     !! mostly we do this to exclude atoms who are involved in short
1191     !! range interactions (bonds, bends, torsions), but we also need
1192     !! to exclude some overcounted interactions that result from
1193     !! the parallel decomposition
1194    
1195 gezelter 1610 #ifdef IS_MPI
1196 gezelter 2204 !! in MPI, we have to look up the unique IDs for each atom
1197     unique_id_1 = AtomRowToGlobal(atom1)
1198 gezelter 1610 #else
1199 gezelter 2204 !! in the normal loop, the atom numbers are unique
1200     unique_id_1 = atom1
1201 gezelter 1610 #endif
1202 gezelter 2204
1203     !! We were called with only one atom, so just check the global exclude
1204     !! list for this atom
1205     if (.not. present(atom2)) then
1206     do i = 1, nExcludes_global
1207     if (excludesGlobal(i) == unique_id_1) then
1208     skip_it = .true.
1209     return
1210     end if
1211     end do
1212     return
1213     end if
1214    
1215 gezelter 1610 #ifdef IS_MPI
1216 gezelter 2204 unique_id_2 = AtomColToGlobal(atom2)
1217 gezelter 1610 #else
1218 gezelter 2204 unique_id_2 = atom2
1219 gezelter 1610 #endif
1220 gezelter 2204
1221 gezelter 1610 #ifdef IS_MPI
1222 gezelter 2204 !! this situation should only arise in MPI simulations
1223     if (unique_id_1 == unique_id_2) then
1224     skip_it = .true.
1225     return
1226     end if
1227    
1228     !! this prevents us from doing the pair on multiple processors
1229     if (unique_id_1 < unique_id_2) then
1230     if (mod(unique_id_1 + unique_id_2,2) == 0) then
1231     skip_it = .true.
1232     return
1233     endif
1234     else
1235     if (mod(unique_id_1 + unique_id_2,2) == 1) then
1236     skip_it = .true.
1237     return
1238     endif
1239     endif
1240 gezelter 1610 #endif
1241 gezelter 2204
1242     !! the rest of these situations can happen in all simulations:
1243     do i = 1, nExcludes_global
1244     if ((excludesGlobal(i) == unique_id_1) .or. &
1245     (excludesGlobal(i) == unique_id_2)) then
1246     skip_it = .true.
1247     return
1248     endif
1249     enddo
1250    
1251     do i = 1, nSkipsForAtom(atom1)
1252     if (skipsForAtom(atom1, i) .eq. unique_id_2) then
1253     skip_it = .true.
1254     return
1255     endif
1256     end do
1257    
1258     return
1259     end function skipThisPair
1260    
1261     function FF_UsesDirectionalAtoms() result(doesit)
1262     logical :: doesit
1263     doesit = FF_uses_DirectionalAtoms .or. FF_uses_Dipoles .or. &
1264     FF_uses_Quadrupoles .or. FF_uses_Sticky .or. &
1265 chrisfen 2220 FF_uses_StickyPower .or. FF_uses_GayBerne .or. FF_uses_Shapes
1266 gezelter 2204 end function FF_UsesDirectionalAtoms
1267    
1268     function FF_RequiresPrepairCalc() result(doesit)
1269     logical :: doesit
1270     doesit = FF_uses_EAM
1271     end function FF_RequiresPrepairCalc
1272    
1273     function FF_RequiresPostpairCalc() result(doesit)
1274     logical :: doesit
1275     doesit = FF_uses_RF
1276     end function FF_RequiresPostpairCalc
1277    
1278 gezelter 1610 #ifdef PROFILE
1279 gezelter 2204 function getforcetime() result(totalforcetime)
1280     real(kind=dp) :: totalforcetime
1281     totalforcetime = forcetime
1282     end function getforcetime
1283 gezelter 1610 #endif
1284    
1285 gezelter 2204 !! This cleans componets of force arrays belonging only to fortran
1286    
1287     subroutine add_stress_tensor(dpair, fpair)
1288    
1289     real( kind = dp ), dimension(3), intent(in) :: dpair, fpair
1290    
1291     ! because the d vector is the rj - ri vector, and
1292     ! because fx, fy, fz are the force on atom i, we need a
1293     ! negative sign here:
1294    
1295     tau_Temp(1) = tau_Temp(1) - dpair(1) * fpair(1)
1296     tau_Temp(2) = tau_Temp(2) - dpair(1) * fpair(2)
1297     tau_Temp(3) = tau_Temp(3) - dpair(1) * fpair(3)
1298     tau_Temp(4) = tau_Temp(4) - dpair(2) * fpair(1)
1299     tau_Temp(5) = tau_Temp(5) - dpair(2) * fpair(2)
1300     tau_Temp(6) = tau_Temp(6) - dpair(2) * fpair(3)
1301     tau_Temp(7) = tau_Temp(7) - dpair(3) * fpair(1)
1302     tau_Temp(8) = tau_Temp(8) - dpair(3) * fpair(2)
1303     tau_Temp(9) = tau_Temp(9) - dpair(3) * fpair(3)
1304    
1305     virial_Temp = virial_Temp + &
1306     (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
1307    
1308     end subroutine add_stress_tensor
1309    
1310 gezelter 1610 end module doForces