ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/doForces.F90
Revision: 939
Committed: Thu Apr 20 18:24:24 2006 UTC (19 years ago) by gezelter
Original Path: trunk/src/UseTheForce/doForces.F90
File size: 52097 byte(s)
Log Message:
Complete rewrite of spline code and everything that uses it.

File Contents

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