ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/doForces.F90
Revision: 699
Committed: Wed Oct 26 23:31:18 2005 UTC (19 years, 6 months ago) by chrisfen
File size: 47220 byte(s)
Log Message:
reaction field looks to be working now

File Contents

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