ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2355
Committed: Wed Oct 12 18:59:16 2005 UTC (18 years, 8 months ago) by chuckv
File size: 46853 byte(s)
Log Message:
Breaky Breaky: Add Support for seperating potential contributions.

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