ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/doForces.F90
Revision: 610
Committed: Sun Sep 18 20:45:38 2005 UTC (19 years, 7 months ago) by chrisfen
Original Path: trunk/src/UseTheForce/doForces.F90
File size: 42022 byte(s)
Log Message:
reaction field seems to work now, still need to do some testing...

File Contents

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