ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/doForces.F90
Revision: 1483
Committed: Tue Jul 27 21:17:31 2010 UTC (14 years, 9 months ago) by gezelter
File size: 67530 byte(s)
Log Message:
Added GB module to the C++ side, got rid of it on the fortran side.

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

Properties

Name Value
svn:keywords Author Id Revision Date