ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2350
Committed: Mon Oct 10 21:20:46 2005 UTC (18 years, 9 months ago) by chuckv
File size: 46244 byte(s)
Log Message:
Fixed MPI bug with Row and Column indexing of groupToGtype. We now have two seperate maps groupToGtypeRow and groupToGypeCol. GroupToGtypeCol is a pointer that is set to groupToGtypeRow in the single processor version.

File Contents

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