ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/doForces.F90
Revision: 1352
Committed: Tue May 26 13:57:29 2009 UTC (15 years, 11 months ago) by gezelter
Original Path: trunk/src/UseTheForce/doForces.F90
File size: 61079 byte(s)
Log Message:
fixed one more bug for charge - non-charge interactions

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