ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/doForces.F90
Revision: 2503
Committed: Thu Dec 8 22:04:40 2005 UTC (18 years, 7 months ago) by gezelter
File size: 49383 byte(s)
Log Message:
added info printout

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