ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2787
Committed: Mon Jun 5 18:24:45 2006 UTC (18 years, 1 month ago) by gezelter
File size: 51331 byte(s)
Log Message:
Massive changes for GB code with multiple ellipsoid types (a la
Cleaver's paper).

Also, changes in hydrodynamics code to make HydroProp a somewhat
smarter class (rather than just a struct).

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