ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/doForces.F90
Revision: 1133
Committed: Tue May 22 19:30:27 2007 UTC (17 years, 11 months ago) by chuckv
Original Path: trunk/src/UseTheForce/doForces.F90
File size: 58118 byte(s)
Log Message:
Fixed pot in parallel.

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