ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2512
Committed: Thu Dec 15 21:43:16 2005 UTC (18 years, 7 months ago) by gezelter
File size: 49421 byte(s)
Log Message:
fixed a cutoff bug

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