ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2533
Committed: Fri Dec 30 23:15:59 2005 UTC (18 years, 6 months ago) by chuckv
File size: 49844 byte(s)
Log Message:
More Sutton-Chen bug fixes.

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