ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/doForces.F90
Revision: 1467
Committed: Sat Jul 17 15:33:03 2010 UTC (14 years, 9 months ago) by gezelter
File size: 66459 byte(s)
Log Message:
well, it compiles, but still segfaults

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date