ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/doForces.F90
Revision: 1386
Committed: Fri Oct 23 18:41:09 2009 UTC (15 years, 6 months ago) by gezelter
File size: 63517 byte(s)
Log Message:
removing MPI responsibilities from the lowest level force routines.  This is
in preparation for migrating these routines (LJ, electrostatic, eam, suttonchen,
gay-berne, sticky) to C++

File Contents

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