ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2375
Committed: Mon Oct 17 19:12:45 2005 UTC (18 years, 9 months ago) by gezelter
File size: 46903 byte(s)
Log Message:
changing GB architecture

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