ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/doForces.F90
Revision: 2267
Committed: Fri Jul 29 17:34:06 2005 UTC (18 years, 11 months ago) by tim
File size: 42096 byte(s)
Log Message:
fix a bug which does not update me_i and me_j correctly

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