ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/doForces.F90
Revision: 1286
Committed: Wed Sep 10 17:57:55 2008 UTC (16 years, 8 months ago) by gezelter
Original Path: trunk/src/UseTheForce/doForces.F90
File size: 59750 byte(s)
Log Message:
Added support for scaled 1-2, 1-3 and 1-4 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     !! 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 gezelter 1286 !! @version $Id: doForces.F90,v 1.95 2008-09-10 17:57:55 gezelter Exp $, $Date: 2008-09-10 17:57:55 $, $Name: not supported by cvs2svn $, $Revision: 1.95 $
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 571 do i = 1, nAtypes
312 gezelter 582 if (SimHasAtype(i)) then
313 gezelter 575 call getElementProperty(atypes, i, "is_LennardJones", i_is_LJ)
314     call getElementProperty(atypes, i, "is_Electrostatic", i_is_Elect)
315     call getElementProperty(atypes, i, "is_Sticky", i_is_Sticky)
316     call getElementProperty(atypes, i, "is_StickyPower", i_is_StickyP)
317     call getElementProperty(atypes, i, "is_GayBerne", i_is_GB)
318     call getElementProperty(atypes, i, "is_EAM", i_is_EAM)
319     call getElementProperty(atypes, i, "is_Shape", i_is_Shape)
320 chuckv 831 call getElementProperty(atypes, i, "is_SC", i_is_SC)
321 chuckv 599
322 chrisfen 618 if (haveDefaultCutoffs) then
323     atypeMaxCutoff(i) = defaultRcut
324     else
325     if (i_is_LJ) then
326     thisRcut = getSigma(i) * 2.5_dp
327     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
328     endif
329     if (i_is_Elect) then
330     thisRcut = defaultRcut
331     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
332     endif
333     if (i_is_Sticky) then
334     thisRcut = getStickyCut(i)
335     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
336     endif
337     if (i_is_StickyP) then
338     thisRcut = getStickyPowerCut(i)
339     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
340     endif
341     if (i_is_GB) then
342     thisRcut = getGayBerneCut(i)
343     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
344     endif
345     if (i_is_EAM) then
346     thisRcut = getEAMCut(i)
347     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
348     endif
349     if (i_is_Shape) then
350     thisRcut = getShapeCut(i)
351     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
352     endif
353 chuckv 831 if (i_is_SC) then
354     thisRcut = getSCCut(i)
355     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
356     endif
357 gezelter 575 endif
358 gezelter 762
359 gezelter 575 if (atypeMaxCutoff(i).gt.biggestAtypeCutoff) then
360     biggestAtypeCutoff = atypeMaxCutoff(i)
361     endif
362 chrisfen 618
363 gezelter 574 endif
364 gezelter 575 enddo
365 gezelter 581
366 gezelter 575 istart = 1
367 chuckv 651 jstart = 1
368 gezelter 575 #ifdef IS_MPI
369     iend = nGroupsInRow
370 chuckv 651 jend = nGroupsInCol
371 gezelter 575 #else
372     iend = nGroups
373 chuckv 651 jend = nGroups
374 gezelter 575 #endif
375 gezelter 582
376 gezelter 581 !! allocate the groupToGtype and gtypeMaxCutoff here.
377 chuckv 651 if(.not.allocated(groupToGtypeRow)) then
378     ! allocate(groupToGtype(iend))
379     allocate(groupToGtypeRow(iend))
380     else
381     deallocate(groupToGtypeRow)
382     allocate(groupToGtypeRow(iend))
383 chuckv 583 endif
384 chuckv 651 if(.not.allocated(groupMaxCutoffRow)) then
385     allocate(groupMaxCutoffRow(iend))
386     else
387     deallocate(groupMaxCutoffRow)
388     allocate(groupMaxCutoffRow(iend))
389     end if
390    
391     if(.not.allocated(gtypeMaxCutoffRow)) then
392     allocate(gtypeMaxCutoffRow(iend))
393     else
394     deallocate(gtypeMaxCutoffRow)
395     allocate(gtypeMaxCutoffRow(iend))
396     endif
397    
398    
399     #ifdef IS_MPI
400     ! We only allocate new storage if we are in MPI because Ncol /= Nrow
401 chuckv 652 if(.not.associated(groupToGtypeCol)) then
402 chuckv 651 allocate(groupToGtypeCol(jend))
403     else
404     deallocate(groupToGtypeCol)
405     allocate(groupToGtypeCol(jend))
406     end if
407    
408 tim 833 if(.not.associated(groupMaxCutoffCol)) then
409     allocate(groupMaxCutoffCol(jend))
410 chuckv 651 else
411 tim 833 deallocate(groupMaxCutoffCol)
412     allocate(groupMaxCutoffCol(jend))
413 chuckv 651 end if
414 chuckv 652 if(.not.associated(gtypeMaxCutoffCol)) then
415 chuckv 651 allocate(gtypeMaxCutoffCol(jend))
416     else
417     deallocate(gtypeMaxCutoffCol)
418     allocate(gtypeMaxCutoffCol(jend))
419     end if
420    
421     groupMaxCutoffCol = 0.0_dp
422     gtypeMaxCutoffCol = 0.0_dp
423    
424     #endif
425     groupMaxCutoffRow = 0.0_dp
426     gtypeMaxCutoffRow = 0.0_dp
427    
428    
429 gezelter 582 !! first we do a single loop over the cutoff groups to find the
430     !! largest cutoff for any atypes present in this group. We also
431     !! create gtypes at this point.
432    
433 gezelter 960 tol = 1.0e-6_dp
434 chuckv 651 nGroupTypesRow = 0
435 tim 833 nGroupTypesCol = 0
436 gezelter 581 do i = istart, iend
437 gezelter 575 n_in_i = groupStartRow(i+1) - groupStartRow(i)
438 chuckv 651 groupMaxCutoffRow(i) = 0.0_dp
439 gezelter 581 do ia = groupStartRow(i), groupStartRow(i+1)-1
440     atom1 = groupListRow(ia)
441 gezelter 575 #ifdef IS_MPI
442 gezelter 581 me_i = atid_row(atom1)
443 gezelter 575 #else
444 gezelter 581 me_i = atid(atom1)
445     #endif
446 chuckv 651 if (atypeMaxCutoff(me_i).gt.groupMaxCutoffRow(i)) then
447     groupMaxCutoffRow(i)=atypeMaxCutoff(me_i)
448 gezelter 587 endif
449 gezelter 581 enddo
450 chuckv 651 if (nGroupTypesRow.eq.0) then
451     nGroupTypesRow = nGroupTypesRow + 1
452     gtypeMaxCutoffRow(nGroupTypesRow) = groupMaxCutoffRow(i)
453     groupToGtypeRow(i) = nGroupTypesRow
454 gezelter 581 else
455 gezelter 587 GtypeFound = .false.
456 chuckv 651 do g = 1, nGroupTypesRow
457     if ( abs(groupMaxCutoffRow(i) - gtypeMaxCutoffRow(g)).lt.tol) then
458     groupToGtypeRow(i) = g
459 gezelter 587 GtypeFound = .true.
460 gezelter 581 endif
461     enddo
462 gezelter 587 if (.not.GtypeFound) then
463 chuckv 651 nGroupTypesRow = nGroupTypesRow + 1
464     gtypeMaxCutoffRow(nGroupTypesRow) = groupMaxCutoffRow(i)
465     groupToGtypeRow(i) = nGroupTypesRow
466 gezelter 587 endif
467 gezelter 581 endif
468 gezelter 587 enddo
469    
470 chuckv 651 #ifdef IS_MPI
471     do j = jstart, jend
472     n_in_j = groupStartCol(j+1) - groupStartCol(j)
473     groupMaxCutoffCol(j) = 0.0_dp
474     do ja = groupStartCol(j), groupStartCol(j+1)-1
475     atom1 = groupListCol(ja)
476    
477     me_j = atid_col(atom1)
478    
479     if (atypeMaxCutoff(me_j).gt.groupMaxCutoffCol(j)) then
480     groupMaxCutoffCol(j)=atypeMaxCutoff(me_j)
481     endif
482     enddo
483    
484     if (nGroupTypesCol.eq.0) then
485     nGroupTypesCol = nGroupTypesCol + 1
486     gtypeMaxCutoffCol(nGroupTypesCol) = groupMaxCutoffCol(j)
487     groupToGtypeCol(j) = nGroupTypesCol
488     else
489     GtypeFound = .false.
490     do g = 1, nGroupTypesCol
491     if ( abs(groupMaxCutoffCol(j) - gtypeMaxCutoffCol(g)).lt.tol) then
492     groupToGtypeCol(j) = g
493     GtypeFound = .true.
494     endif
495     enddo
496     if (.not.GtypeFound) then
497     nGroupTypesCol = nGroupTypesCol + 1
498     gtypeMaxCutoffCol(nGroupTypesCol) = groupMaxCutoffCol(j)
499     groupToGtypeCol(j) = nGroupTypesCol
500     endif
501     endif
502     enddo
503    
504     #else
505     ! Set pointers to information we just found
506     nGroupTypesCol = nGroupTypesRow
507     groupToGtypeCol => groupToGtypeRow
508     gtypeMaxCutoffCol => gtypeMaxCutoffRow
509     groupMaxCutoffCol => groupMaxCutoffRow
510     #endif
511    
512 gezelter 581 !! allocate the gtypeCutoffMap here.
513 chuckv 651 allocate(gtypeCutoffMap(nGroupTypesRow,nGroupTypesCol))
514 gezelter 581 !! then we do a double loop over all the group TYPES to find the cutoff
515     !! map between groups of two types
516 chuckv 651 tradRcut = max(maxval(gtypeMaxCutoffRow),maxval(gtypeMaxCutoffCol))
517    
518 gezelter 762 do i = 1, nGroupTypesRow
519 chuckv 651 do j = 1, nGroupTypesCol
520 gezelter 576
521 gezelter 581 select case(cutoffPolicy)
522 gezelter 582 case(TRADITIONAL_CUTOFF_POLICY)
523 chuckv 651 thisRcut = tradRcut
524 gezelter 582 case(MIX_CUTOFF_POLICY)
525 chuckv 651 thisRcut = 0.5_dp * (gtypeMaxCutoffRow(i) + gtypeMaxCutoffCol(j))
526 gezelter 582 case(MAX_CUTOFF_POLICY)
527 chuckv 651 thisRcut = max(gtypeMaxCutoffRow(i), gtypeMaxCutoffCol(j))
528 gezelter 582 case default
529     call handleError("createGtypeCutoffMap", "Unknown Cutoff Policy")
530     return
531     end select
532     gtypeCutoffMap(i,j)%rcut = thisRcut
533 gezelter 762
534     if (thisRcut.gt.largestRcut) largestRcut = thisRcut
535    
536 gezelter 582 gtypeCutoffMap(i,j)%rcutsq = thisRcut*thisRcut
537 gezelter 585
538 gezelter 762 if (.not.haveSkinThickness) then
539     skinThickness = 1.0_dp
540     endif
541    
542     gtypeCutoffMap(i,j)%rlistsq = (thisRcut + skinThickness)**2
543    
544 chrisfen 618 ! sanity check
545    
546     if (haveDefaultCutoffs) then
547     if (abs(gtypeCutoffMap(i,j)%rcut - defaultRcut).gt.0.0001) then
548     call handleError("createGtypeCutoffMap", "user-specified rCut does not match computed group Cutoff")
549     endif
550     endif
551 gezelter 581 enddo
552     enddo
553 gezelter 762
554 chuckv 651 if(allocated(gtypeMaxCutoffRow)) deallocate(gtypeMaxCutoffRow)
555     if(allocated(groupMaxCutoffRow)) deallocate(groupMaxCutoffRow)
556     if(allocated(atypeMaxCutoff)) deallocate(atypeMaxCutoff)
557     #ifdef IS_MPI
558     if(associated(groupMaxCutoffCol)) deallocate(groupMaxCutoffCol)
559     if(associated(gtypeMaxCutoffCol)) deallocate(gtypeMaxCutoffCol)
560     #endif
561     groupMaxCutoffCol => null()
562     gtypeMaxCutoffCol => null()
563    
564 gezelter 581 haveGtypeCutoffMap = .true.
565 chrisfen 596 end subroutine createGtypeCutoffMap
566 chrisfen 578
567 chrisfen 1129 subroutine setCutoffs(defRcut, defRsw, defSP, defSF)
568 chrisfen 596
569 gezelter 762 real(kind=dp),intent(in) :: defRcut, defRsw
570 chrisfen 1129 logical, intent(in) :: defSP, defSF
571 gezelter 762 character(len = statusMsgSize) :: errMsg
572     integer :: localError
573    
574 chrisfen 596 defaultRcut = defRcut
575     defaultRsw = defRsw
576 gezelter 762
577 chrisfen 1129 defaultDoShiftPot = defSP
578     defaultDoShiftFrc = defSF
579    
580 gezelter 762 if (abs(defaultRcut-defaultRsw) .lt. 0.0001) then
581 chrisfen 1129 if (defaultDoShiftFrc) then
582     write(errMsg, *) &
583     'cutoffRadius and switchingRadius are set to the', newline &
584     // tab, 'same value. OOPSE will use shifted force', newline &
585     // tab, 'potentials instead of switching functions.'
586    
587     call handleInfo("setCutoffs", errMsg)
588     else
589     write(errMsg, *) &
590     'cutoffRadius and switchingRadius are set to the', newline &
591     // tab, 'same value. OOPSE will use shifted', newline &
592     // tab, 'potentials instead of switching functions.'
593    
594     call handleInfo("setCutoffs", errMsg)
595    
596     defaultDoShiftPot = .true.
597     endif
598    
599 gezelter 762 endif
600 gezelter 939
601 gezelter 762 localError = 0
602 chrisfen 1129 call setLJDefaultCutoff( defaultRcut, defaultDoShiftPot, &
603     defaultDoShiftFrc )
604 gezelter 813 call setElectrostaticCutoffRadius( defaultRcut, defaultRsw )
605 gezelter 938 call setCutoffEAM( defaultRcut )
606     call setCutoffSC( defaultRcut )
607 chuckv 1162 call setMnMDefaultCutoff( defaultRcut, defaultDoShiftPot, &
608     defaultDoShiftFrc )
609 gezelter 939 call set_switch(defaultRsw, defaultRcut)
610 gezelter 889 call setHmatDangerousRcutValue(defaultRcut)
611 gezelter 939
612 chrisfen 618 haveDefaultCutoffs = .true.
613 gezelter 813 haveGtypeCutoffMap = .false.
614 gezelter 939
615 gezelter 762 end subroutine setCutoffs
616 chrisfen 596
617 gezelter 762 subroutine cWasLame()
618    
619     VisitCutoffsAfterComputing = .true.
620     return
621    
622     end subroutine cWasLame
623    
624 chrisfen 596 subroutine setCutoffPolicy(cutPolicy)
625 gezelter 762
626 chrisfen 596 integer, intent(in) :: cutPolicy
627 gezelter 762
628 chrisfen 596 cutoffPolicy = cutPolicy
629 gezelter 762 haveCutoffPolicy = .true.
630 gezelter 813 haveGtypeCutoffMap = .false.
631 gezelter 762
632 gezelter 576 end subroutine setCutoffPolicy
633 gezelter 1126
634 chrisfen 998 subroutine setBoxDipole()
635    
636     do_box_dipole = .true.
637    
638     end subroutine setBoxDipole
639    
640     subroutine getBoxDipole( box_dipole )
641    
642     real(kind=dp), intent(inout), dimension(3) :: box_dipole
643    
644     box_dipole = boxDipole
645    
646     end subroutine getBoxDipole
647    
648 gezelter 762 subroutine setElectrostaticMethod( thisESM )
649    
650     integer, intent(in) :: thisESM
651    
652     electrostaticSummationMethod = thisESM
653     haveElectrostaticSummationMethod = .true.
654 gezelter 574
655 gezelter 762 end subroutine setElectrostaticMethod
656    
657     subroutine setSkinThickness( thisSkin )
658 gezelter 574
659 gezelter 762 real(kind=dp), intent(in) :: thisSkin
660    
661     skinThickness = thisSkin
662 gezelter 813 haveSkinThickness = .true.
663     haveGtypeCutoffMap = .false.
664 gezelter 762
665     end subroutine setSkinThickness
666    
667     subroutine setSimVariables()
668     SIM_uses_DirectionalAtoms = SimUsesDirectionalAtoms()
669     SIM_uses_EAM = SimUsesEAM()
670     SIM_requires_postpair_calc = SimRequiresPostpairCalc()
671     SIM_requires_prepair_calc = SimRequiresPrepairCalc()
672     SIM_uses_PBC = SimUsesPBC()
673 chuckv 841 SIM_uses_SC = SimUsesSC()
674 gezelter 1126 SIM_uses_AtomicVirial = SimUsesAtomicVirial()
675 chrisfen 998
676 gezelter 762 haveSIMvariables = .true.
677    
678     return
679     end subroutine setSimVariables
680 gezelter 117
681     subroutine doReadyCheck(error)
682     integer, intent(out) :: error
683     integer :: myStatus
684    
685     error = 0
686 chrisfen 532
687 gezelter 571 if (.not. haveInteractionHash) then
688 gezelter 762 call createInteractionHash()
689 gezelter 117 endif
690    
691 gezelter 571 if (.not. haveGtypeCutoffMap) then
692 gezelter 762 call createGtypeCutoffMap()
693 gezelter 571 endif
694    
695 gezelter 762 if (VisitCutoffsAfterComputing) then
696 gezelter 939 call set_switch(largestRcut, largestRcut)
697 gezelter 889 call setHmatDangerousRcutValue(largestRcut)
698 gezelter 938 call setCutoffEAM(largestRcut)
699     call setCutoffSC(largestRcut)
700     VisitCutoffsAfterComputing = .false.
701 gezelter 762 endif
702    
703 gezelter 117 if (.not. haveSIMvariables) then
704     call setSimVariables()
705     endif
706    
707     if (.not. haveNeighborList) then
708     write(default_error, *) 'neighbor list has not been initialized in doForces!'
709     error = -1
710     return
711     end if
712 gezelter 939
713 gezelter 117 if (.not. haveSaneForceField) then
714     write(default_error, *) 'Force Field is not sane in doForces!'
715     error = -1
716     return
717     end if
718 gezelter 939
719 gezelter 117 #ifdef IS_MPI
720     if (.not. isMPISimSet()) then
721     write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
722     error = -1
723     return
724     endif
725     #endif
726     return
727     end subroutine doReadyCheck
728    
729 chrisfen 532
730 gezelter 762 subroutine init_FF(thisStat)
731 gezelter 117
732     integer, intent(out) :: thisStat
733     integer :: my_status, nMatches
734     integer, pointer :: MatchList(:) => null()
735    
736     !! assume things are copacetic, unless they aren't
737     thisStat = 0
738    
739     !! init_FF is called *after* all of the atom types have been
740     !! defined in atype_module using the new_atype subroutine.
741     !!
742     !! this will scan through the known atypes and figure out what
743     !! interactions are used by the force field.
744 chrisfen 532
745 gezelter 141 FF_uses_DirectionalAtoms = .false.
746     FF_uses_Dipoles = .false.
747     FF_uses_GayBerne = .false.
748 gezelter 117 FF_uses_EAM = .false.
749 chuckv 834 FF_uses_SC = .false.
750 chrisfen 532
751 gezelter 141 call getMatchingElementList(atypes, "is_Directional", .true., &
752     nMatches, MatchList)
753     if (nMatches .gt. 0) FF_uses_DirectionalAtoms = .true.
754    
755     call getMatchingElementList(atypes, "is_Dipole", .true., &
756     nMatches, MatchList)
757 gezelter 571 if (nMatches .gt. 0) FF_uses_Dipoles = .true.
758 chrisfen 523
759 gezelter 141 call getMatchingElementList(atypes, "is_GayBerne", .true., &
760     nMatches, MatchList)
761 gezelter 571 if (nMatches .gt. 0) FF_uses_GayBerne = .true.
762 chrisfen 532
763 gezelter 117 call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList)
764     if (nMatches .gt. 0) FF_uses_EAM = .true.
765 chrisfen 532
766 chuckv 834 call getMatchingElementList(atypes, "is_SC", .true., nMatches, MatchList)
767     if (nMatches .gt. 0) FF_uses_SC = .true.
768 gezelter 141
769 chuckv 834
770 gezelter 117 haveSaneForceField = .true.
771 chrisfen 532
772 gezelter 117 if (FF_uses_EAM) then
773 chrisfen 532 call init_EAM_FF(my_status)
774 gezelter 117 if (my_status /= 0) then
775     write(default_error, *) "init_EAM_FF returned a bad status"
776     thisStat = -1
777     haveSaneForceField = .false.
778     return
779     end if
780     endif
781    
782     if (.not. haveNeighborList) then
783     !! Create neighbor lists
784     call expandNeighborList(nLocal, my_status)
785     if (my_Status /= 0) then
786     write(default_error,*) "SimSetup: ExpandNeighborList returned error."
787     thisStat = -1
788     return
789     endif
790     haveNeighborList = .true.
791 chrisfen 532 endif
792    
793 gezelter 117 end subroutine init_FF
794    
795 chrisfen 532
796 gezelter 117 !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
797     !------------------------------------------------------------->
798 gezelter 1285 subroutine do_force_loop(q, q_group, A, eFrame, f, t, tau, pot, particle_pot, &
799 gezelter 117 do_pot_c, do_stress_c, error)
800     !! Position array provided by C, dimensioned by getNlocal
801     real ( kind = dp ), dimension(3, nLocal) :: q
802     !! molecular center-of-mass position array
803     real ( kind = dp ), dimension(3, nGroups) :: q_group
804     !! Rotation Matrix for each long range particle in simulation.
805     real( kind = dp), dimension(9, nLocal) :: A
806     !! Unit vectors for dipoles (lab frame)
807 gezelter 246 real( kind = dp ), dimension(9,nLocal) :: eFrame
808 gezelter 117 !! Force array provided by C, dimensioned by getNlocal
809     real ( kind = dp ), dimension(3,nLocal) :: f
810     !! Torsion array provided by C, dimensioned by getNlocal
811     real( kind = dp ), dimension(3,nLocal) :: t
812    
813     !! Stress Tensor
814     real( kind = dp), dimension(9) :: tau
815 gezelter 662 real ( kind = dp ),dimension(LR_POT_TYPES) :: pot
816 chuckv 1245 real( kind = dp ), dimension(nLocal) :: particle_pot
817 gezelter 117 logical ( kind = 2) :: do_pot_c, do_stress_c
818     logical :: do_pot
819     logical :: do_stress
820     logical :: in_switching_region
821     #ifdef IS_MPI
822 gezelter 662 real( kind = DP ), dimension(LR_POT_TYPES) :: pot_local
823 gezelter 117 integer :: nAtomsInRow
824     integer :: nAtomsInCol
825     integer :: nprocs
826     integer :: nGroupsInRow
827     integer :: nGroupsInCol
828     #endif
829     integer :: natoms
830     logical :: update_nlist
831     integer :: i, j, jstart, jend, jnab
832     integer :: istart, iend
833     integer :: ia, jb, atom1, atom2
834     integer :: nlist
835 gezelter 1126 real( kind = DP ) :: ratmsq, rgrpsq, rgrp, rag, vpair, vij
836 gezelter 117 real( kind = DP ) :: sw, dswdr, swderiv, mf
837 chrisfen 699 real( kind = DP ) :: rVal
838 gezelter 1126 real(kind=dp),dimension(3) :: d_atm, d_grp, fpair, fij, fg, dag
839     real(kind=dp) :: rfpot, mu_i
840 gezelter 762 real(kind=dp):: rCut
841 gezelter 117 integer :: me_i, me_j, n_in_i, n_in_j
842     logical :: is_dp_i
843     integer :: neighborListSize
844     integer :: listerror, error
845     integer :: localError
846     integer :: propPack_i, propPack_j
847     integer :: loopStart, loopEnd, loop
848 gezelter 571 integer :: iHash
849 gezelter 1286 integer :: i1, topoDist
850 chrisfen 532
851 chrisfen 998 !! the variables for the box dipole moment
852     #ifdef IS_MPI
853     integer :: pChgCount_local
854     integer :: nChgCount_local
855     real(kind=dp) :: pChg_local
856     real(kind=dp) :: nChg_local
857     real(kind=dp), dimension(3) :: pChgPos_local
858     real(kind=dp), dimension(3) :: nChgPos_local
859     real(kind=dp), dimension(3) :: dipVec_local
860     #endif
861     integer :: pChgCount
862     integer :: nChgCount
863     real(kind=dp) :: pChg
864     real(kind=dp) :: nChg
865     real(kind=dp) :: chg_value
866     real(kind=dp), dimension(3) :: pChgPos
867     real(kind=dp), dimension(3) :: nChgPos
868     real(kind=dp), dimension(3) :: dipVec
869     real(kind=dp), dimension(3) :: chgVec
870    
871     !! initialize box dipole variables
872     if (do_box_dipole) then
873     #ifdef IS_MPI
874     pChg_local = 0.0_dp
875     nChg_local = 0.0_dp
876     pChgCount_local = 0
877     nChgCount_local = 0
878     do i=1, 3
879     pChgPos_local = 0.0_dp
880     nChgPos_local = 0.0_dp
881     dipVec_local = 0.0_dp
882     enddo
883     #endif
884     pChg = 0.0_dp
885     nChg = 0.0_dp
886     pChgCount = 0
887     nChgCount = 0
888     chg_value = 0.0_dp
889    
890     do i=1, 3
891     pChgPos(i) = 0.0_dp
892     nChgPos(i) = 0.0_dp
893     dipVec(i) = 0.0_dp
894     chgVec(i) = 0.0_dp
895     boxDipole(i) = 0.0_dp
896     enddo
897     endif
898    
899 gezelter 117 !! initialize local variables
900 chrisfen 532
901 gezelter 117 #ifdef IS_MPI
902     pot_local = 0.0_dp
903     nAtomsInRow = getNatomsInRow(plan_atom_row)
904     nAtomsInCol = getNatomsInCol(plan_atom_col)
905     nGroupsInRow = getNgroupsInRow(plan_group_row)
906     nGroupsInCol = getNgroupsInCol(plan_group_col)
907     #else
908     natoms = nlocal
909     #endif
910 chrisfen 532
911 gezelter 117 call doReadyCheck(localError)
912     if ( localError .ne. 0 ) then
913     call handleError("do_force_loop", "Not Initialized")
914     error = -1
915     return
916     end if
917     call zero_work_arrays()
918 chrisfen 532
919 gezelter 117 do_pot = do_pot_c
920     do_stress = do_stress_c
921 chrisfen 532
922 gezelter 117 ! Gather all information needed by all force loops:
923 chrisfen 532
924 gezelter 117 #ifdef IS_MPI
925 chrisfen 532
926 gezelter 117 call gather(q, q_Row, plan_atom_row_3d)
927     call gather(q, q_Col, plan_atom_col_3d)
928    
929     call gather(q_group, q_group_Row, plan_group_row_3d)
930     call gather(q_group, q_group_Col, plan_group_col_3d)
931 chrisfen 532
932 gezelter 141 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
933 gezelter 246 call gather(eFrame, eFrame_Row, plan_atom_row_rotation)
934     call gather(eFrame, eFrame_Col, plan_atom_col_rotation)
935 chrisfen 532
936 gezelter 117 call gather(A, A_Row, plan_atom_row_rotation)
937     call gather(A, A_Col, plan_atom_col_rotation)
938     endif
939 chrisfen 532
940 gezelter 117 #endif
941 chrisfen 532
942 gezelter 117 !! Begin force loop timing:
943     #ifdef PROFILE
944     call cpu_time(forceTimeInitial)
945     nloops = nloops + 1
946     #endif
947 chrisfen 532
948 gezelter 117 loopEnd = PAIR_LOOP
949     if (FF_RequiresPrepairCalc() .and. SIM_requires_prepair_calc) then
950     loopStart = PREPAIR_LOOP
951     else
952     loopStart = PAIR_LOOP
953     endif
954    
955     do loop = loopStart, loopEnd
956    
957     ! See if we need to update neighbor lists
958     ! (but only on the first time through):
959     if (loop .eq. loopStart) then
960     #ifdef IS_MPI
961 gezelter 762 call checkNeighborList(nGroupsInRow, q_group_row, skinThickness, &
962 chrisfen 532 update_nlist)
963 gezelter 117 #else
964 gezelter 762 call checkNeighborList(nGroups, q_group, skinThickness, &
965 chrisfen 532 update_nlist)
966 gezelter 117 #endif
967     endif
968 chrisfen 532
969 gezelter 117 if (update_nlist) then
970     !! save current configuration and construct neighbor list
971     #ifdef IS_MPI
972     call saveNeighborList(nGroupsInRow, q_group_row)
973     #else
974     call saveNeighborList(nGroups, q_group)
975     #endif
976     neighborListSize = size(list)
977     nlist = 0
978     endif
979 chrisfen 532
980 gezelter 117 istart = 1
981     #ifdef IS_MPI
982     iend = nGroupsInRow
983     #else
984     iend = nGroups - 1
985     #endif
986     outer: do i = istart, iend
987    
988     if (update_nlist) point(i) = nlist + 1
989 chrisfen 532
990 gezelter 117 n_in_i = groupStartRow(i+1) - groupStartRow(i)
991 chrisfen 532
992 gezelter 117 if (update_nlist) then
993     #ifdef IS_MPI
994     jstart = 1
995     jend = nGroupsInCol
996     #else
997     jstart = i+1
998     jend = nGroups
999     #endif
1000     else
1001     jstart = point(i)
1002     jend = point(i+1) - 1
1003     ! make sure group i has neighbors
1004     if (jstart .gt. jend) cycle outer
1005     endif
1006 chrisfen 532
1007 gezelter 117 do jnab = jstart, jend
1008     if (update_nlist) then
1009     j = jnab
1010     else
1011     j = list(jnab)
1012     endif
1013    
1014     #ifdef IS_MPI
1015 chuckv 567 me_j = atid_col(j)
1016 gezelter 117 call get_interatomic_vector(q_group_Row(:,i), &
1017     q_group_Col(:,j), d_grp, rgrpsq)
1018     #else
1019 chuckv 567 me_j = atid(j)
1020 gezelter 117 call get_interatomic_vector(q_group(:,i), &
1021     q_group(:,j), d_grp, rgrpsq)
1022 chrisfen 618 #endif
1023 gezelter 117
1024 chuckv 651 if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rListsq) then
1025 gezelter 117 if (update_nlist) then
1026     nlist = nlist + 1
1027 chrisfen 532
1028 gezelter 117 if (nlist > neighborListSize) then
1029     #ifdef IS_MPI
1030     call expandNeighborList(nGroupsInRow, listerror)
1031     #else
1032     call expandNeighborList(nGroups, listerror)
1033     #endif
1034     if (listerror /= 0) then
1035     error = -1
1036     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
1037     return
1038     end if
1039     neighborListSize = size(list)
1040     endif
1041 chrisfen 532
1042 gezelter 117 list(nlist) = j
1043     endif
1044 gezelter 939
1045 chrisfen 708 if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCutsq) then
1046 chrisfen 532
1047 gezelter 762 rCut = gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCut
1048 chrisfen 708 if (loop .eq. PAIR_LOOP) then
1049 gezelter 960 vij = 0.0_dp
1050 gezelter 938 fij(1) = 0.0_dp
1051     fij(2) = 0.0_dp
1052     fij(3) = 0.0_dp
1053 chrisfen 708 endif
1054    
1055 gezelter 939 call get_switch(rgrpsq, sw, dswdr,rgrp, in_switching_region)
1056 chrisfen 708
1057     n_in_j = groupStartCol(j+1) - groupStartCol(j)
1058    
1059     do ia = groupStartRow(i), groupStartRow(i+1)-1
1060 chrisfen 703
1061 chrisfen 708 atom1 = groupListRow(ia)
1062    
1063     inner: do jb = groupStartCol(j), groupStartCol(j+1)-1
1064    
1065     atom2 = groupListCol(jb)
1066    
1067     if (skipThisPair(atom1, atom2)) cycle inner
1068    
1069     if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
1070 gezelter 938 d_atm(1) = d_grp(1)
1071     d_atm(2) = d_grp(2)
1072     d_atm(3) = d_grp(3)
1073 chrisfen 708 ratmsq = rgrpsq
1074     else
1075 gezelter 117 #ifdef IS_MPI
1076 chrisfen 708 call get_interatomic_vector(q_Row(:,atom1), &
1077     q_Col(:,atom2), d_atm, ratmsq)
1078 gezelter 117 #else
1079 chrisfen 708 call get_interatomic_vector(q(:,atom1), &
1080     q(:,atom2), d_atm, ratmsq)
1081 gezelter 117 #endif
1082 gezelter 1286 endif
1083    
1084     topoDist = getTopoDistance(atom1, atom2)
1085    
1086 chrisfen 708 if (loop .eq. PREPAIR_LOOP) then
1087 gezelter 117 #ifdef IS_MPI
1088 chrisfen 708 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
1089 gezelter 762 rgrpsq, d_grp, rCut, do_pot, do_stress, &
1090 chrisfen 708 eFrame, A, f, t, pot_local)
1091 gezelter 117 #else
1092 chrisfen 708 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
1093 gezelter 762 rgrpsq, d_grp, rCut, do_pot, do_stress, &
1094 chrisfen 708 eFrame, A, f, t, pot)
1095 gezelter 117 #endif
1096 chrisfen 708 else
1097 gezelter 117 #ifdef IS_MPI
1098 chrisfen 708 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
1099     do_pot, eFrame, A, f, t, pot_local, vpair, &
1100 gezelter 1286 fpair, d_grp, rgrp, rCut, topoDist)
1101 chuckv 1245 ! particle_pot will be accumulated from row & column
1102     ! arrays later
1103 gezelter 117 #else
1104 chrisfen 708 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
1105 chuckv 1245 do_pot, eFrame, A, f, t, pot, vpair, &
1106 gezelter 1286 fpair, d_grp, rgrp, rCut, topoDist)
1107 chuckv 1245 particle_pot(atom1) = particle_pot(atom1) + vpair*sw
1108     particle_pot(atom2) = particle_pot(atom2) + vpair*sw
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    
1296    
1297 gezelter 117 endif
1298     #endif
1299 chrisfen 532
1300 chrisfen 691 if (SIM_requires_postpair_calc) then
1301 chrisfen 695 do i = 1, nlocal
1302    
1303     ! we loop only over the local atoms, so we don't need row and column
1304     ! lookups for the types
1305 chrisfen 699
1306 chrisfen 691 me_i = atid(i)
1307    
1308 chrisfen 695 ! is the atom electrostatic? See if it would have an
1309     ! electrostatic interaction with itself
1310     iHash = InteractionHash(me_i,me_i)
1311 chrisfen 699
1312 chrisfen 691 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1313 gezelter 117 #ifdef IS_MPI
1314 chrisfen 703 call self_self(i, eFrame, pot_local(ELECTROSTATIC_POT), &
1315 chrisfen 695 t, do_pot)
1316 gezelter 117 #else
1317 chrisfen 703 call self_self(i, eFrame, pot(ELECTROSTATIC_POT), &
1318 chrisfen 695 t, do_pot)
1319 gezelter 117 #endif
1320 chrisfen 691 endif
1321 chrisfen 699
1322 chrisfen 703
1323 chrisfen 708 if (electrostaticSummationMethod.eq.REACTION_FIELD) then
1324 chrisfen 699
1325 chrisfen 703 ! loop over the excludes to accumulate RF stuff we've
1326     ! left out of the normal pair loop
1327    
1328     do i1 = 1, nSkipsForAtom(i)
1329     j = skipsForAtom(i, i1)
1330    
1331     ! prevent overcounting of the skips
1332     if (i.lt.j) then
1333 gezelter 939 call get_interatomic_vector(q(:,i), q(:,j), d_atm, ratmsq)
1334 gezelter 960 rVal = sqrt(ratmsq)
1335 gezelter 939 call get_switch(ratmsq, sw, dswdr, rVal,in_switching_region)
1336 chrisfen 699 #ifdef IS_MPI
1337 gezelter 1286 call rf_self_excludes(i, j, sw, 1.0_dp, eFrame, d_atm, rVal, &
1338 chrisfen 703 vpair, pot_local(ELECTROSTATIC_POT), f, t, do_pot)
1339 chrisfen 699 #else
1340 gezelter 1286 call rf_self_excludes(i, j, sw, 1.0_dp, eFrame, d_atm, rVal, &
1341 chrisfen 703 vpair, pot(ELECTROSTATIC_POT), f, t, do_pot)
1342 chrisfen 699 #endif
1343 chrisfen 703 endif
1344     enddo
1345 chrisfen 708 endif
1346 chrisfen 998
1347     if (do_box_dipole) then
1348     #ifdef IS_MPI
1349     call accumulate_box_dipole(i, eFrame, q(:,i), pChg_local, &
1350     nChg_local, pChgPos_local, nChgPos_local, dipVec_local, &
1351     pChgCount_local, nChgCount_local)
1352     #else
1353     call accumulate_box_dipole(i, eFrame, q(:,i), pChg, nChg, &
1354     pChgPos, nChgPos, dipVec, pChgCount, nChgCount)
1355     #endif
1356     endif
1357 chrisfen 703 enddo
1358 gezelter 117 endif
1359 chrisfen 998
1360 gezelter 117 #ifdef IS_MPI
1361     if (do_pot) then
1362 gezelter 962 #ifdef SINGLE_PRECISION
1363     call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_real,mpi_sum, &
1364     mpi_comm_world,mpi_err)
1365     #else
1366 chrisfen 998 call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_double_precision, &
1367     mpi_sum, mpi_comm_world,mpi_err)
1368 gezelter 962 #endif
1369 gezelter 117 endif
1370 gezelter 1126
1371 chrisfen 998 if (do_box_dipole) then
1372    
1373     #ifdef SINGLE_PRECISION
1374     call mpi_allreduce(pChg_local, pChg, 1, mpi_real, mpi_sum, &
1375     mpi_comm_world, mpi_err)
1376     call mpi_allreduce(nChg_local, nChg, 1, mpi_real, mpi_sum, &
1377     mpi_comm_world, mpi_err)
1378     call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer, mpi_sum,&
1379     mpi_comm_world, mpi_err)
1380     call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer, mpi_sum,&
1381     mpi_comm_world, mpi_err)
1382     call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_real, mpi_sum, &
1383     mpi_comm_world, mpi_err)
1384     call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_real, mpi_sum, &
1385     mpi_comm_world, mpi_err)
1386     call mpi_allreduce(dipVec_local, dipVec, 3, mpi_real, mpi_sum, &
1387     mpi_comm_world, mpi_err)
1388 gezelter 117 #else
1389 chrisfen 998 call mpi_allreduce(pChg_local, pChg, 1, mpi_double_precision, mpi_sum, &
1390     mpi_comm_world, mpi_err)
1391     call mpi_allreduce(nChg_local, nChg, 1, mpi_double_precision, mpi_sum, &
1392     mpi_comm_world, mpi_err)
1393     call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer,&
1394     mpi_sum, mpi_comm_world, mpi_err)
1395     call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer,&
1396     mpi_sum, mpi_comm_world, mpi_err)
1397     call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_double_precision, &
1398     mpi_sum, mpi_comm_world, mpi_err)
1399     call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_double_precision, &
1400     mpi_sum, mpi_comm_world, mpi_err)
1401     call mpi_allreduce(dipVec_local, dipVec, 3, mpi_double_precision, &
1402     mpi_sum, mpi_comm_world, mpi_err)
1403     #endif
1404    
1405     endif
1406 chrisfen 695
1407 gezelter 117 #endif
1408 chrisfen 998
1409     if (do_box_dipole) then
1410     ! first load the accumulated dipole moment (if dipoles were present)
1411     boxDipole(1) = dipVec(1)
1412     boxDipole(2) = dipVec(2)
1413     boxDipole(3) = dipVec(3)
1414    
1415     ! now include the dipole moment due to charges
1416     ! use the lesser of the positive and negative charge totals
1417     if (nChg .le. pChg) then
1418     chg_value = nChg
1419     else
1420     chg_value = pChg
1421     endif
1422    
1423     ! find the average positions
1424     if (pChgCount .gt. 0 .and. nChgCount .gt. 0) then
1425     pChgPos = pChgPos / pChgCount
1426     nChgPos = nChgPos / nChgCount
1427     endif
1428    
1429     ! dipole is from the negative to the positive (physics notation)
1430     chgVec(1) = pChgPos(1) - nChgPos(1)
1431     chgVec(2) = pChgPos(2) - nChgPos(2)
1432     chgVec(3) = pChgPos(3) - nChgPos(3)
1433    
1434     boxDipole(1) = boxDipole(1) + chgVec(1) * chg_value
1435     boxDipole(2) = boxDipole(2) + chgVec(2) * chg_value
1436     boxDipole(3) = boxDipole(3) + chgVec(3) * chg_value
1437    
1438     endif
1439    
1440 gezelter 117 end subroutine do_force_loop
1441 chrisfen 532
1442 gezelter 117 subroutine do_pair(i, j, rijsq, d, sw, do_pot, &
1443 gezelter 1286 eFrame, A, f, t, pot, vpair, fpair, d_grp, r_grp, rCut, topoDist)
1444 gezelter 117
1445 chuckv 656 real( kind = dp ) :: vpair, sw
1446 gezelter 662 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1447 chuckv 1245 real( kind = dp ), dimension(nLocal) :: particle_pot
1448 gezelter 117 real( kind = dp ), dimension(3) :: fpair
1449     real( kind = dp ), dimension(nLocal) :: mfact
1450 gezelter 246 real( kind = dp ), dimension(9,nLocal) :: eFrame
1451 gezelter 117 real( kind = dp ), dimension(9,nLocal) :: A
1452     real( kind = dp ), dimension(3,nLocal) :: f
1453     real( kind = dp ), dimension(3,nLocal) :: t
1454    
1455     logical, intent(inout) :: do_pot
1456     integer, intent(in) :: i, j
1457     real ( kind = dp ), intent(inout) :: rijsq
1458 chrisfen 695 real ( kind = dp ), intent(inout) :: r_grp
1459 gezelter 117 real ( kind = dp ), intent(inout) :: d(3)
1460 chrisfen 695 real ( kind = dp ), intent(inout) :: d_grp(3)
1461 gezelter 762 real ( kind = dp ), intent(inout) :: rCut
1462 gezelter 1286 integer, intent(inout) :: topoDist
1463     real ( kind = dp ) :: r, pair_pot, vdwMult, electroMult
1464 gezelter 939 real ( kind = dp ) :: a_k, b_k, c_k, d_k, dx
1465 gezelter 117 integer :: me_i, me_j
1466 gezelter 939 integer :: k
1467 gezelter 117
1468 gezelter 571 integer :: iHash
1469 gezelter 560
1470 chrisfen 942 r = sqrt(rijsq)
1471    
1472 gezelter 960 vpair = 0.0_dp
1473     fpair(1:3) = 0.0_dp
1474 gezelter 117
1475     #ifdef IS_MPI
1476     me_i = atid_row(i)
1477     me_j = atid_col(j)
1478     #else
1479     me_i = atid(i)
1480     me_j = atid(j)
1481     #endif
1482 gezelter 202
1483 gezelter 571 iHash = InteractionHash(me_i, me_j)
1484 gezelter 1286 vdwMult = vdwScale(topoDist)
1485     electroMult = electrostaticScale(topoDist)
1486 chrisfen 703
1487     if ( iand(iHash, LJ_PAIR).ne.0 ) then
1488 gezelter 1286 call do_lj_pair(i, j, d, r, rijsq, rcut, sw, vdwMult, vpair, fpair, &
1489 chrisfen 703 pot(VDW_POT), f, do_pot)
1490 gezelter 117 endif
1491 chrisfen 532
1492 chrisfen 703 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1493 gezelter 1286 call doElectrostaticPair(i, j, d, r, rijsq, rcut, sw, electroMult, &
1494     vpair, fpair, pot(ELECTROSTATIC_POT), eFrame, f, t, do_pot)
1495 chrisfen 703 endif
1496    
1497     if ( iand(iHash, STICKY_PAIR).ne.0 ) then
1498     call do_sticky_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1499     pot(HB_POT), A, f, t, do_pot)
1500     endif
1501    
1502     if ( iand(iHash, STICKYPOWER_PAIR).ne.0 ) then
1503     call do_sticky_power_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1504     pot(HB_POT), A, f, t, do_pot)
1505     endif
1506    
1507     if ( iand(iHash, GAYBERNE_PAIR).ne.0 ) then
1508 gezelter 1286 call do_gb_pair(i, j, d, r, rijsq, sw, vdwMult, vpair, fpair, &
1509 chrisfen 703 pot(VDW_POT), A, f, t, do_pot)
1510     endif
1511    
1512     if ( iand(iHash, GAYBERNE_LJ).ne.0 ) then
1513 gezelter 1286 call do_gb_pair(i, j, d, r, rijsq, sw, vdwMult, vpair, fpair, &
1514 chrisfen 703 pot(VDW_POT), A, f, t, do_pot)
1515     endif
1516    
1517     if ( iand(iHash, EAM_PAIR).ne.0 ) then
1518     call do_eam_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1519     pot(METALLIC_POT), f, do_pot)
1520     endif
1521    
1522     if ( iand(iHash, SHAPE_PAIR).ne.0 ) then
1523     call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1524     pot(VDW_POT), A, f, t, do_pot)
1525     endif
1526    
1527     if ( iand(iHash, SHAPE_LJ).ne.0 ) then
1528     call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1529     pot(VDW_POT), A, f, t, do_pot)
1530     endif
1531 chuckv 733
1532     if ( iand(iHash, SC_PAIR).ne.0 ) then
1533 gezelter 762 call do_SC_pair(i, j, d, r, rijsq, rcut, sw, vpair, fpair, &
1534 chuckv 733 pot(METALLIC_POT), f, do_pot)
1535     endif
1536 chrisfen 703
1537 gezelter 1174 if ( iand(iHash, MNM_PAIR).ne.0 ) then
1538 gezelter 1286 call do_mnm_pair(i, j, d, r, rijsq, rcut, sw, vdwMult, vpair, fpair, &
1539 gezelter 1174 pot(VDW_POT), A, f, t, do_pot)
1540     endif
1541    
1542 gezelter 1286
1543     if (topoDist .ne. 0) then
1544     write(*,*) i, j, vpair, fpair(1), fpair(2), fpair(3)
1545     endif
1546    
1547 gezelter 117 end subroutine do_pair
1548    
1549 gezelter 762 subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, rCut, &
1550 gezelter 246 do_pot, do_stress, eFrame, A, f, t, pot)
1551 gezelter 117
1552 chuckv 656 real( kind = dp ) :: sw
1553 gezelter 662 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1554 chrisfen 532 real( kind = dp ), dimension(9,nLocal) :: eFrame
1555     real (kind=dp), dimension(9,nLocal) :: A
1556     real (kind=dp), dimension(3,nLocal) :: f
1557     real (kind=dp), dimension(3,nLocal) :: t
1558 gezelter 117
1559 chrisfen 532 logical, intent(inout) :: do_pot, do_stress
1560     integer, intent(in) :: i, j
1561 gezelter 762 real ( kind = dp ), intent(inout) :: rijsq, rcijsq, rCut
1562 chrisfen 532 real ( kind = dp ) :: r, rc
1563     real ( kind = dp ), intent(inout) :: d(3), dc(3)
1564    
1565 gezelter 571 integer :: me_i, me_j, iHash
1566 chrisfen 532
1567 chrisfen 942 r = sqrt(rijsq)
1568    
1569 gezelter 117 #ifdef IS_MPI
1570 chrisfen 532 me_i = atid_row(i)
1571     me_j = atid_col(j)
1572 gezelter 117 #else
1573 chrisfen 532 me_i = atid(i)
1574     me_j = atid(j)
1575 gezelter 117 #endif
1576 chrisfen 532
1577 gezelter 571 iHash = InteractionHash(me_i, me_j)
1578 chrisfen 532
1579 gezelter 571 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1580 gezelter 762 call calc_EAM_prepair_rho(i, j, d, r, rijsq)
1581 chrisfen 532 endif
1582 chuckv 733
1583     if ( iand(iHash, SC_PAIR).ne.0 ) then
1584 gezelter 762 call calc_SC_prepair_rho(i, j, d, r, rijsq, rcut )
1585 chuckv 733 endif
1586 gezelter 560
1587 chrisfen 532 end subroutine do_prepair
1588    
1589    
1590 gezelter 1285 subroutine do_preforce(nlocal, pot, particle_pot)
1591 chrisfen 532 integer :: nlocal
1592 gezelter 662 real( kind = dp ),dimension(LR_POT_TYPES) :: pot
1593 gezelter 1285 real( kind = dp ),dimension(nlocal) :: particle_pot
1594 chrisfen 532
1595     if (FF_uses_EAM .and. SIM_uses_EAM) then
1596 gezelter 1285 call calc_EAM_preforce_Frho(nlocal, pot(METALLIC_POT), particle_pot)
1597 chrisfen 532 endif
1598 chuckv 733 if (FF_uses_SC .and. SIM_uses_SC) then
1599 gezelter 1285 call calc_SC_preforce_Frho(nlocal, pot(METALLIC_POT), particle_pot)
1600 chuckv 733 endif
1601 chrisfen 532 end subroutine do_preforce
1602    
1603    
1604     subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
1605    
1606     real (kind = dp), dimension(3) :: q_i
1607     real (kind = dp), dimension(3) :: q_j
1608     real ( kind = dp ), intent(out) :: r_sq
1609     real( kind = dp ) :: d(3), scaled(3)
1610     integer i
1611    
1612 gezelter 938 d(1) = q_j(1) - q_i(1)
1613     d(2) = q_j(2) - q_i(2)
1614     d(3) = q_j(3) - q_i(3)
1615 chrisfen 532
1616     ! Wrap back into periodic box if necessary
1617     if ( SIM_uses_PBC ) then
1618    
1619     if( .not.boxIsOrthorhombic ) then
1620     ! calc the scaled coordinates.
1621 gezelter 939 ! scaled = matmul(HmatInv, d)
1622 chrisfen 532
1623 gezelter 939 scaled(1) = HmatInv(1,1)*d(1) + HmatInv(1,2)*d(2) + HmatInv(1,3)*d(3)
1624     scaled(2) = HmatInv(2,1)*d(1) + HmatInv(2,2)*d(2) + HmatInv(2,3)*d(3)
1625     scaled(3) = HmatInv(3,1)*d(1) + HmatInv(3,2)*d(2) + HmatInv(3,3)*d(3)
1626    
1627 chrisfen 532 ! wrap the scaled coordinates
1628    
1629 gezelter 960 scaled(1) = scaled(1) - anint(scaled(1), kind=dp)
1630     scaled(2) = scaled(2) - anint(scaled(2), kind=dp)
1631     scaled(3) = scaled(3) - anint(scaled(3), kind=dp)
1632 chrisfen 532
1633     ! calc the wrapped real coordinates from the wrapped scaled
1634     ! coordinates
1635 gezelter 939 ! d = matmul(Hmat,scaled)
1636     d(1)= Hmat(1,1)*scaled(1) + Hmat(1,2)*scaled(2) + Hmat(1,3)*scaled(3)
1637     d(2)= Hmat(2,1)*scaled(1) + Hmat(2,2)*scaled(2) + Hmat(2,3)*scaled(3)
1638     d(3)= Hmat(3,1)*scaled(1) + Hmat(3,2)*scaled(2) + Hmat(3,3)*scaled(3)
1639 chrisfen 532
1640     else
1641     ! calc the scaled coordinates.
1642    
1643 gezelter 938 scaled(1) = d(1) * HmatInv(1,1)
1644     scaled(2) = d(2) * HmatInv(2,2)
1645     scaled(3) = d(3) * HmatInv(3,3)
1646    
1647     ! wrap the scaled coordinates
1648    
1649 gezelter 960 scaled(1) = scaled(1) - anint(scaled(1), kind=dp)
1650     scaled(2) = scaled(2) - anint(scaled(2), kind=dp)
1651     scaled(3) = scaled(3) - anint(scaled(3), kind=dp)
1652 chrisfen 532
1653 gezelter 938 ! calc the wrapped real coordinates from the wrapped scaled
1654     ! coordinates
1655 chrisfen 532
1656 gezelter 938 d(1) = scaled(1)*Hmat(1,1)
1657     d(2) = scaled(2)*Hmat(2,2)
1658     d(3) = scaled(3)*Hmat(3,3)
1659 chrisfen 532
1660     endif
1661    
1662     endif
1663    
1664 gezelter 938 r_sq = d(1)*d(1) + d(2)*d(2) + d(3)*d(3)
1665 chrisfen 532
1666     end subroutine get_interatomic_vector
1667    
1668     subroutine zero_work_arrays()
1669    
1670 gezelter 117 #ifdef IS_MPI
1671    
1672 chrisfen 532 q_Row = 0.0_dp
1673     q_Col = 0.0_dp
1674    
1675     q_group_Row = 0.0_dp
1676     q_group_Col = 0.0_dp
1677    
1678     eFrame_Row = 0.0_dp
1679     eFrame_Col = 0.0_dp
1680    
1681     A_Row = 0.0_dp
1682     A_Col = 0.0_dp
1683    
1684     f_Row = 0.0_dp
1685     f_Col = 0.0_dp
1686     f_Temp = 0.0_dp
1687    
1688     t_Row = 0.0_dp
1689     t_Col = 0.0_dp
1690     t_Temp = 0.0_dp
1691    
1692     pot_Row = 0.0_dp
1693     pot_Col = 0.0_dp
1694     pot_Temp = 0.0_dp
1695    
1696 gezelter 117 #endif
1697 chrisfen 532
1698     if (FF_uses_EAM .and. SIM_uses_EAM) then
1699     call clean_EAM()
1700     endif
1701    
1702     end subroutine zero_work_arrays
1703    
1704     function skipThisPair(atom1, atom2) result(skip_it)
1705     integer, intent(in) :: atom1
1706     integer, intent(in), optional :: atom2
1707     logical :: skip_it
1708     integer :: unique_id_1, unique_id_2
1709     integer :: me_i,me_j
1710     integer :: i
1711    
1712     skip_it = .false.
1713    
1714     !! there are a number of reasons to skip a pair or a particle
1715     !! mostly we do this to exclude atoms who are involved in short
1716     !! range interactions (bonds, bends, torsions), but we also need
1717     !! to exclude some overcounted interactions that result from
1718     !! the parallel decomposition
1719    
1720 gezelter 117 #ifdef IS_MPI
1721 chrisfen 532 !! in MPI, we have to look up the unique IDs for each atom
1722     unique_id_1 = AtomRowToGlobal(atom1)
1723     unique_id_2 = AtomColToGlobal(atom2)
1724     !! this situation should only arise in MPI simulations
1725     if (unique_id_1 == unique_id_2) then
1726     skip_it = .true.
1727     return
1728     end if
1729    
1730     !! this prevents us from doing the pair on multiple processors
1731     if (unique_id_1 < unique_id_2) then
1732     if (mod(unique_id_1 + unique_id_2,2) == 0) then
1733     skip_it = .true.
1734     return
1735     endif
1736     else
1737     if (mod(unique_id_1 + unique_id_2,2) == 1) then
1738     skip_it = .true.
1739     return
1740     endif
1741     endif
1742 gezelter 1286 #else
1743     !! in the normal loop, the atom numbers are unique
1744     unique_id_1 = atom1
1745     unique_id_2 = atom2
1746 gezelter 117 #endif
1747 gezelter 1286
1748 chrisfen 532 do i = 1, nSkipsForAtom(atom1)
1749     if (skipsForAtom(atom1, i) .eq. unique_id_2) then
1750     skip_it = .true.
1751     return
1752     endif
1753     end do
1754    
1755     return
1756     end function skipThisPair
1757    
1758 gezelter 1286 function getTopoDistance(atom1, atom2) result(topoDist)
1759     integer, intent(in) :: atom1
1760     integer, intent(in) :: atom2
1761     integer :: topoDist
1762     integer :: unique_id_2
1763     integer :: i
1764    
1765     #ifdef IS_MPI
1766     unique_id_2 = AtomColToGlobal(atom2)
1767     #else
1768     unique_id_2 = atom2
1769     #endif
1770    
1771     ! zero is default for unconnected (i.e. normal) pair interactions
1772    
1773     topoDist = 0
1774    
1775     do i = 1, nTopoPairsForAtom(atom1)
1776     if (toposForAtom(atom1, i) .eq. unique_id_2) then
1777     topoDist = topoDistance(atom1, i)
1778     return
1779     endif
1780     end do
1781    
1782     return
1783     end function getTopoDistance
1784    
1785 chrisfen 532 function FF_UsesDirectionalAtoms() result(doesit)
1786     logical :: doesit
1787 gezelter 571 doesit = FF_uses_DirectionalAtoms
1788 chrisfen 532 end function FF_UsesDirectionalAtoms
1789    
1790     function FF_RequiresPrepairCalc() result(doesit)
1791     logical :: doesit
1792 chuckv 1162 doesit = FF_uses_EAM .or. FF_uses_SC
1793 chrisfen 532 end function FF_RequiresPrepairCalc
1794    
1795 gezelter 117 #ifdef PROFILE
1796 chrisfen 532 function getforcetime() result(totalforcetime)
1797     real(kind=dp) :: totalforcetime
1798     totalforcetime = forcetime
1799     end function getforcetime
1800 gezelter 117 #endif
1801    
1802 chrisfen 532 !! This cleans componets of force arrays belonging only to fortran
1803    
1804 gezelter 1126 subroutine add_stress_tensor(dpair, fpair, tau)
1805 chrisfen 532
1806     real( kind = dp ), dimension(3), intent(in) :: dpair, fpair
1807 gezelter 1126 real( kind = dp ), dimension(9), intent(inout) :: tau
1808 chrisfen 532
1809     ! because the d vector is the rj - ri vector, and
1810     ! because fx, fy, fz are the force on atom i, we need a
1811     ! negative sign here:
1812    
1813 gezelter 1126 tau(1) = tau(1) - dpair(1) * fpair(1)
1814     tau(2) = tau(2) - dpair(1) * fpair(2)
1815     tau(3) = tau(3) - dpair(1) * fpair(3)
1816     tau(4) = tau(4) - dpair(2) * fpair(1)
1817     tau(5) = tau(5) - dpair(2) * fpair(2)
1818     tau(6) = tau(6) - dpair(2) * fpair(3)
1819     tau(7) = tau(7) - dpair(3) * fpair(1)
1820     tau(8) = tau(8) - dpair(3) * fpair(2)
1821     tau(9) = tau(9) - dpair(3) * fpair(3)
1822 chrisfen 532
1823     end subroutine add_stress_tensor
1824    
1825 gezelter 117 end module doForces