ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/doForces.F90
Revision: 2390
Committed: Wed Oct 19 19:24:40 2005 UTC (18 years, 8 months ago) by chrisfen
File size: 46069 byte(s)
Log Message:
Still had some globals toUpper problems - these changes should fix those...

File Contents

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