ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/doForces.F90
Revision: 1419
Committed: Fri Mar 26 17:25:54 2010 UTC (15 years, 1 month ago) by gezelter
File size: 67634 byte(s)
Log Message:
Fixed a typo in EAM, cleaned up EAM and sutton chen a bit

File Contents

# Content
1 !!
2 !! Copyright (c) 2005, 2009 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. Redistributions of source code must retain the above copyright
10 !! notice, this list of conditions and the following disclaimer.
11 !!
12 !! 2. Redistributions in binary form must reproduce the above copyright
13 !! notice, this list of conditions and the following disclaimer in the
14 !! documentation and/or other materials provided with the
15 !! distribution.
16 !!
17 !! This software is provided "AS IS," without a warranty of any
18 !! kind. All express or implied conditions, representations and
19 !! warranties, including any implied warranty of merchantability,
20 !! fitness for a particular purpose or non-infringement, are hereby
21 !! excluded. The University of Notre Dame and its licensors shall not
22 !! be liable for any damages suffered by licensee as a result of
23 !! using, modifying or distributing the software or its
24 !! derivatives. In no event will the University of Notre Dame or its
25 !! licensors be liable for any lost revenue, profit or data, or for
26 !! direct, indirect, special, consequential, incidental or punitive
27 !! damages, however caused and regardless of the theory of liability,
28 !! arising out of the use of or inability to use software, even if the
29 !! University of Notre Dame has been advised of the possibility of
30 !! such damages.
31 !!
32 !! SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33 !! research, please cite the appropriate papers when you publish your
34 !! work. Good starting points are:
35 !!
36 !! [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37 !! [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38 !! [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39 !! [4] Vardeman & Gezelter, in progress (2009).
40 !!
41
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.106 2009-11-25 20:01:57 gezelter Exp $, $Date: 2009-11-25 20:01:57 $, $Name: not supported by cvs2svn $, $Revision: 1.106 $
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. OpenMD 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. OpenMD 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
782 if (.not. haveNeighborList) then
783 !! Create neighbor lists
784 call expandNeighborList(nLocal, my_status)
785 if (my_Status /= 0) then
786 write(default_error,*) "SimSetup: ExpandNeighborList returned error."
787 thisStat = -1
788 return
789 endif
790 haveNeighborList = .true.
791 endif
792
793 end subroutine init_FF
794
795
796 !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
797 !------------------------------------------------------------->
798 subroutine do_force_loop(q, q_group, A, eFrame, f, t, tau, pot, particle_pot, &
799 do_pot_c, do_stress_c, error)
800 !! Position array provided by C, dimensioned by getNlocal
801 real ( kind = dp ), dimension(3, nLocal) :: q
802 !! molecular center-of-mass position array
803 real ( kind = dp ), dimension(3, nGroups) :: q_group
804 !! Rotation Matrix for each long range particle in simulation.
805 real( kind = dp), dimension(9, nLocal) :: A
806 !! Unit vectors for dipoles (lab frame)
807 real( kind = dp ), dimension(9,nLocal) :: eFrame
808 !! Force array provided by C, dimensioned by getNlocal
809 real ( kind = dp ), dimension(3,nLocal) :: f
810 !! Torsion array provided by C, dimensioned by getNlocal
811 real( kind = dp ), dimension(3,nLocal) :: t
812
813 !! Stress Tensor
814 real( kind = dp), dimension(9) :: tau
815 real ( kind = dp ),dimension(LR_POT_TYPES) :: pot
816 real( kind = dp ), dimension(nLocal) :: particle_pot
817 logical ( kind = 2) :: do_pot_c, do_stress_c
818 logical :: do_pot
819 logical :: do_stress
820 logical :: in_switching_region
821 #ifdef IS_MPI
822 real( kind = DP ), dimension(LR_POT_TYPES) :: pot_local
823 integer :: nAtomsInRow
824 integer :: nAtomsInCol
825 integer :: nprocs
826 integer :: nGroupsInRow
827 integer :: nGroupsInCol
828 #endif
829 integer :: natoms
830 logical :: update_nlist
831 integer :: i, j, jstart, jend, jnab
832 integer :: istart, iend
833 integer :: ia, jb, atom1, atom2
834 integer :: nlist
835 real( kind = DP ) :: ratmsq, rgrpsq, rgrp, rag, vpair, vij
836 real( kind = DP ) :: sw, dswdr, swderiv, mf
837 real( kind = DP ) :: rVal
838 real(kind=dp),dimension(3) :: d_atm, d_grp, fpair, fij, fg, dag
839 real(kind=dp) :: rfpot, mu_i
840 real(kind=dp):: rCut
841 integer :: me_i, me_j, n_in_i, n_in_j, iG, j1
842 logical :: is_dp_i
843 integer :: neighborListSize
844 integer :: listerror, error
845 integer :: localError
846 integer :: propPack_i, propPack_j
847 integer :: loopStart, loopEnd, loop
848 integer :: iHash, jHash
849 integer :: i1, topoDist
850
851 !! the variables for the box dipole moment
852 #ifdef IS_MPI
853 integer :: pChgCount_local
854 integer :: nChgCount_local
855 real(kind=dp) :: pChg_local
856 real(kind=dp) :: nChg_local
857 real(kind=dp), dimension(3) :: pChgPos_local
858 real(kind=dp), dimension(3) :: nChgPos_local
859 real(kind=dp), dimension(3) :: dipVec_local
860 #endif
861 integer :: pChgCount
862 integer :: nChgCount
863 real(kind=dp) :: pChg
864 real(kind=dp) :: nChg
865 real(kind=dp) :: chg_value
866 real(kind=dp), dimension(3) :: pChgPos
867 real(kind=dp), dimension(3) :: nChgPos
868 real(kind=dp), dimension(3) :: dipVec
869 real(kind=dp), dimension(3) :: chgVec
870 real(kind=dp) :: skch
871
872 !! initialize box dipole variables
873 if (do_box_dipole) then
874 #ifdef IS_MPI
875 pChg_local = 0.0_dp
876 nChg_local = 0.0_dp
877 pChgCount_local = 0
878 nChgCount_local = 0
879 do i=1, 3
880 pChgPos_local = 0.0_dp
881 nChgPos_local = 0.0_dp
882 dipVec_local = 0.0_dp
883 enddo
884 #endif
885 pChg = 0.0_dp
886 nChg = 0.0_dp
887 pChgCount = 0
888 nChgCount = 0
889 chg_value = 0.0_dp
890
891 do i=1, 3
892 pChgPos(i) = 0.0_dp
893 nChgPos(i) = 0.0_dp
894 dipVec(i) = 0.0_dp
895 chgVec(i) = 0.0_dp
896 boxDipole(i) = 0.0_dp
897 enddo
898 endif
899
900 !! initialize local variables
901
902 #ifdef IS_MPI
903 pot_local = 0.0_dp
904 nAtomsInRow = getNatomsInRow(plan_atom_row)
905 nAtomsInCol = getNatomsInCol(plan_atom_col)
906 nGroupsInRow = getNgroupsInRow(plan_group_row)
907 nGroupsInCol = getNgroupsInCol(plan_group_col)
908 #else
909 natoms = nlocal
910 #endif
911
912 call doReadyCheck(localError)
913 if ( localError .ne. 0 ) then
914 call handleError("do_force_loop", "Not Initialized")
915 error = -1
916 return
917 end if
918 call zero_work_arrays()
919
920 do_pot = do_pot_c
921 do_stress = do_stress_c
922
923 ! Gather all information needed by all force loops:
924
925 #ifdef IS_MPI
926
927 call gather(q, q_Row, plan_atom_row_3d)
928 call gather(q, q_Col, plan_atom_col_3d)
929
930 call gather(q_group, q_group_Row, plan_group_row_3d)
931 call gather(q_group, q_group_Col, plan_group_col_3d)
932
933 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
934 call gather(eFrame, eFrame_Row, plan_atom_row_rotation)
935 call gather(eFrame, eFrame_Col, plan_atom_col_rotation)
936
937 call gather(A, A_Row, plan_atom_row_rotation)
938 call gather(A, A_Col, plan_atom_col_rotation)
939 endif
940
941 #endif
942
943 !! Begin force loop timing:
944 #ifdef PROFILE
945 call cpu_time(forceTimeInitial)
946 nloops = nloops + 1
947 #endif
948
949 loopEnd = PAIR_LOOP
950 if (FF_RequiresPrepairCalc() .and. SIM_requires_prepair_calc) then
951 loopStart = PREPAIR_LOOP
952 else
953 loopStart = PAIR_LOOP
954 endif
955
956 do loop = loopStart, loopEnd
957
958 ! See if we need to update neighbor lists
959 ! (but only on the first time through):
960 if (loop .eq. loopStart) then
961 #ifdef IS_MPI
962 call checkNeighborList(nGroupsInRow, q_group_row, skinThickness, &
963 update_nlist)
964 #else
965 call checkNeighborList(nGroups, q_group, skinThickness, &
966 update_nlist)
967 #endif
968 endif
969
970 if (update_nlist) then
971 !! save current configuration and construct neighbor list
972 #ifdef IS_MPI
973 call saveNeighborList(nGroupsInRow, q_group_row)
974 #else
975 call saveNeighborList(nGroups, q_group)
976 #endif
977 neighborListSize = size(list)
978 nlist = 0
979 endif
980
981 istart = 1
982 #ifdef IS_MPI
983 iend = nGroupsInRow
984 #else
985 iend = nGroups - 1
986 #endif
987 outer: do i = istart, iend
988
989 if (update_nlist) point(i) = nlist + 1
990
991 n_in_i = groupStartRow(i+1) - groupStartRow(i)
992
993 if (update_nlist) then
994 #ifdef IS_MPI
995 jstart = 1
996 jend = nGroupsInCol
997 #else
998 jstart = i+1
999 jend = nGroups
1000 #endif
1001 else
1002 jstart = point(i)
1003 jend = point(i+1) - 1
1004 ! make sure group i has neighbors
1005 if (jstart .gt. jend) cycle outer
1006 endif
1007
1008 do jnab = jstart, jend
1009 if (update_nlist) then
1010 j = jnab
1011 else
1012 j = list(jnab)
1013 endif
1014
1015 #ifdef IS_MPI
1016 me_j = atid_col(j)
1017 call get_interatomic_vector(q_group_Row(:,i), &
1018 q_group_Col(:,j), d_grp, rgrpsq)
1019 #else
1020 me_j = atid(j)
1021 call get_interatomic_vector(q_group(:,i), &
1022 q_group(:,j), d_grp, rgrpsq)
1023 #endif
1024
1025 if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rListsq) then
1026 if (update_nlist) then
1027 nlist = nlist + 1
1028
1029 if (nlist > neighborListSize) then
1030 #ifdef IS_MPI
1031 call expandNeighborList(nGroupsInRow, listerror)
1032 #else
1033 call expandNeighborList(nGroups, listerror)
1034 #endif
1035 if (listerror /= 0) then
1036 error = -1
1037 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
1038 return
1039 end if
1040 neighborListSize = size(list)
1041 endif
1042
1043 list(nlist) = j
1044 endif
1045
1046 if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCutsq) then
1047
1048 rCut = gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCut
1049 if (loop .eq. PAIR_LOOP) then
1050 vij = 0.0_dp
1051 fij(1) = 0.0_dp
1052 fij(2) = 0.0_dp
1053 fij(3) = 0.0_dp
1054 endif
1055
1056 call get_switch(rgrpsq, sw, dswdr,rgrp, in_switching_region)
1057
1058 n_in_j = groupStartCol(j+1) - groupStartCol(j)
1059
1060 do ia = groupStartRow(i), groupStartRow(i+1)-1
1061
1062 atom1 = groupListRow(ia)
1063
1064 inner: do jb = groupStartCol(j), groupStartCol(j+1)-1
1065
1066 atom2 = groupListCol(jb)
1067
1068 if (skipThisPair(atom1, atom2)) cycle inner
1069
1070 if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
1071 d_atm(1) = d_grp(1)
1072 d_atm(2) = d_grp(2)
1073 d_atm(3) = d_grp(3)
1074 ratmsq = rgrpsq
1075 else
1076 #ifdef IS_MPI
1077 call get_interatomic_vector(q_Row(:,atom1), &
1078 q_Col(:,atom2), d_atm, ratmsq)
1079 #else
1080 call get_interatomic_vector(q(:,atom1), &
1081 q(:,atom2), d_atm, ratmsq)
1082 #endif
1083 endif
1084
1085 topoDist = getTopoDistance(atom1, atom2)
1086
1087 if (loop .eq. PREPAIR_LOOP) then
1088 #ifdef IS_MPI
1089 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
1090 rgrpsq, d_grp, rCut, do_pot, do_stress, &
1091 eFrame, A, f, t, pot_local)
1092 #else
1093 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
1094 rgrpsq, d_grp, rCut, do_pot, do_stress, &
1095 eFrame, A, f, t, pot)
1096 #endif
1097 else
1098 #ifdef IS_MPI
1099 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
1100 do_pot, eFrame, A, f, t, pot_local, particle_pot, vpair, &
1101 fpair, d_grp, rgrp, rCut, topoDist)
1102 ! particle_pot will be accumulated from row & column
1103 ! arrays later
1104 #else
1105 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
1106 do_pot, eFrame, A, f, t, pot, particle_pot, vpair, &
1107 fpair, d_grp, rgrp, rCut, topoDist)
1108 #endif
1109 vij = vij + vpair
1110 fij(1) = fij(1) + fpair(1)
1111 fij(2) = fij(2) + fpair(2)
1112 fij(3) = fij(3) + fpair(3)
1113 if (do_stress) then
1114 call add_stress_tensor(d_atm, fpair, tau)
1115 endif
1116 endif
1117 enddo inner
1118 enddo
1119
1120 if (loop .eq. PAIR_LOOP) then
1121 if (in_switching_region) then
1122 swderiv = vij*dswdr/rgrp
1123 fg = swderiv*d_grp
1124
1125 fij(1) = fij(1) + fg(1)
1126 fij(2) = fij(2) + fg(2)
1127 fij(3) = fij(3) + fg(3)
1128
1129 if (do_stress .and. (n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
1130 call add_stress_tensor(d_atm, fg, tau)
1131 endif
1132
1133 do ia=groupStartRow(i), groupStartRow(i+1)-1
1134 atom1=groupListRow(ia)
1135 mf = mfactRow(atom1)
1136 ! fg is the force on atom ia due to cutoff group's
1137 ! presence in switching region
1138 fg = swderiv*d_grp*mf
1139 #ifdef IS_MPI
1140 f_Row(1,atom1) = f_Row(1,atom1) + fg(1)
1141 f_Row(2,atom1) = f_Row(2,atom1) + fg(2)
1142 f_Row(3,atom1) = f_Row(3,atom1) + fg(3)
1143 #else
1144 f(1,atom1) = f(1,atom1) + fg(1)
1145 f(2,atom1) = f(2,atom1) + fg(2)
1146 f(3,atom1) = f(3,atom1) + fg(3)
1147 #endif
1148 if (n_in_i .gt. 1) then
1149 if (do_stress.and.SIM_uses_AtomicVirial) then
1150 ! find the distance between the atom and the center of
1151 ! the cutoff group:
1152 #ifdef IS_MPI
1153 call get_interatomic_vector(q_Row(:,atom1), &
1154 q_group_Row(:,i), dag, rag)
1155 #else
1156 call get_interatomic_vector(q(:,atom1), &
1157 q_group(:,i), dag, rag)
1158 #endif
1159 call add_stress_tensor(dag,fg,tau)
1160 endif
1161 endif
1162 enddo
1163
1164 do jb=groupStartCol(j), groupStartCol(j+1)-1
1165 atom2=groupListCol(jb)
1166 mf = mfactCol(atom2)
1167 ! fg is the force on atom jb due to cutoff group's
1168 ! presence in switching region
1169 fg = -swderiv*d_grp*mf
1170 #ifdef IS_MPI
1171 f_Col(1,atom2) = f_Col(1,atom2) + fg(1)
1172 f_Col(2,atom2) = f_Col(2,atom2) + fg(2)
1173 f_Col(3,atom2) = f_Col(3,atom2) + fg(3)
1174 #else
1175 f(1,atom2) = f(1,atom2) + fg(1)
1176 f(2,atom2) = f(2,atom2) + fg(2)
1177 f(3,atom2) = f(3,atom2) + fg(3)
1178 #endif
1179 if (n_in_j .gt. 1) then
1180 if (do_stress.and.SIM_uses_AtomicVirial) then
1181 ! find the distance between the atom and the center of
1182 ! the cutoff group:
1183 #ifdef IS_MPI
1184 call get_interatomic_vector(q_Col(:,atom2), &
1185 q_group_Col(:,j), dag, rag)
1186 #else
1187 call get_interatomic_vector(q(:,atom2), &
1188 q_group(:,j), dag, rag)
1189 #endif
1190 call add_stress_tensor(dag,fg,tau)
1191 endif
1192 endif
1193 enddo
1194 endif
1195 !if (do_stress.and.(.not.SIM_uses_AtomicVirial)) then
1196 ! call add_stress_tensor(d_grp, fij, tau)
1197 !endif
1198 endif
1199 endif
1200 endif
1201 enddo
1202
1203 enddo outer
1204
1205 if (update_nlist) then
1206 #ifdef IS_MPI
1207 point(nGroupsInRow + 1) = nlist + 1
1208 #else
1209 point(nGroups) = nlist + 1
1210 #endif
1211 if (loop .eq. PREPAIR_LOOP) then
1212 ! we just did the neighbor list update on the first
1213 ! pass, so we don't need to do it
1214 ! again on the second pass
1215 update_nlist = .false.
1216 endif
1217 endif
1218
1219 if (loop .eq. PREPAIR_LOOP) then
1220 #ifdef IS_MPI
1221 call do_preforce(nlocal, pot_local, particle_pot)
1222 #else
1223 call do_preforce(nlocal, pot, particle_pot)
1224 #endif
1225 endif
1226
1227 enddo
1228
1229 !! Do timing
1230 #ifdef PROFILE
1231 call cpu_time(forceTimeFinal)
1232 forceTime = forceTime + forceTimeFinal - forceTimeInitial
1233 #endif
1234
1235 #ifdef IS_MPI
1236 !!distribute forces
1237
1238 f_temp = 0.0_dp
1239 call scatter(f_Row,f_temp,plan_atom_row_3d)
1240 do i = 1,nlocal
1241 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1242 end do
1243
1244 f_temp = 0.0_dp
1245 call scatter(f_Col,f_temp,plan_atom_col_3d)
1246 do i = 1,nlocal
1247 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1248 end do
1249
1250 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
1251 t_temp = 0.0_dp
1252 call scatter(t_Row,t_temp,plan_atom_row_3d)
1253 do i = 1,nlocal
1254 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1255 end do
1256 t_temp = 0.0_dp
1257 call scatter(t_Col,t_temp,plan_atom_col_3d)
1258
1259 do i = 1,nlocal
1260 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1261 end do
1262 endif
1263
1264 if (do_pot) then
1265 ! scatter/gather pot_row into the members of my column
1266 do i = 1,LR_POT_TYPES
1267 call scatter(pot_Row(i,:), pot_Temp(i,:), plan_atom_row)
1268 end do
1269 ! scatter/gather pot_local into all other procs
1270 ! add resultant to get total pot
1271 do i = 1, nlocal
1272 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES) &
1273 + pot_Temp(1:LR_POT_TYPES,i)
1274 enddo
1275
1276 do i = 1,LR_POT_TYPES
1277 particle_pot(1:nlocal) = particle_pot(1:nlocal) + pot_Temp(i,1:nlocal)
1278 enddo
1279
1280 pot_Temp = 0.0_DP
1281
1282 do i = 1,LR_POT_TYPES
1283 call scatter(pot_Col(i,:), pot_Temp(i,:), plan_atom_col)
1284 end do
1285
1286 do i = 1, nlocal
1287 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES)&
1288 + pot_Temp(1:LR_POT_TYPES,i)
1289 enddo
1290
1291 do i = 1,LR_POT_TYPES
1292 particle_pot(1:nlocal) = particle_pot(1:nlocal) + pot_Temp(i,1:nlocal)
1293 enddo
1294
1295 ppot_Temp = 0.0_DP
1296
1297 call scatter(ppot_Row(:), ppot_Temp(:), plan_atom_row)
1298 do i = 1, nlocal
1299 particle_pot(i) = particle_pot(i) + ppot_Temp(i)
1300 enddo
1301
1302 ppot_Temp = 0.0_DP
1303
1304 call scatter(ppot_Col(:), ppot_Temp(:), plan_atom_col)
1305 do i = 1, nlocal
1306 particle_pot(i) = particle_pot(i) + ppot_Temp(i)
1307 enddo
1308
1309
1310 endif
1311 #endif
1312
1313 if (SIM_requires_postpair_calc) then
1314 do i = 1, nlocal
1315
1316 ! we loop only over the local atoms, so we don't need row and column
1317 ! lookups for the types
1318
1319 me_i = atid(i)
1320
1321 ! is the atom electrostatic? See if it would have an
1322 ! electrostatic interaction with itself
1323 iHash = InteractionHash(me_i,me_i)
1324
1325 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1326
1327 ! loop over the excludes to accumulate charge in the
1328 ! cutoff sphere that we've left out of the normal pair loop
1329 skch = 0.0_dp
1330
1331 do i1 = 1, nSkipsForLocalAtom(i)
1332 j = skipsForLocalAtom(i, i1)
1333 me_j = atid(j)
1334 jHash = InteractionHash(me_i,me_j)
1335 if ( iand(jHash, ELECTROSTATIC_PAIR).ne.0 ) then
1336 skch = skch + getCharge(me_j)
1337 endif
1338 enddo
1339
1340 #ifdef IS_MPI
1341 call self_self(i, eFrame, skch, pot_local(ELECTROSTATIC_POT), &
1342 t, do_pot)
1343 #else
1344 call self_self(i, eFrame, skch, pot(ELECTROSTATIC_POT), &
1345 t, do_pot)
1346 #endif
1347 endif
1348
1349
1350 if (electrostaticSummationMethod.eq.REACTION_FIELD) then
1351
1352 ! loop over the excludes to accumulate RF stuff we've
1353 ! left out of the normal pair loop
1354
1355 do i1 = 1, nSkipsForLocalAtom(i)
1356 j = skipsForLocalAtom(i, i1)
1357
1358 ! prevent overcounting of the skips
1359 if (i.lt.j) then
1360 call get_interatomic_vector(q(:,i), q(:,j), d_atm, ratmsq)
1361 rVal = sqrt(ratmsq)
1362 call get_switch(ratmsq, sw, dswdr, rVal,in_switching_region)
1363 #ifdef IS_MPI
1364 call rf_self_excludes(i, j, sw, 1.0_dp, eFrame, d_atm, rVal, &
1365 vpair, pot_local(ELECTROSTATIC_POT), f, t, do_pot)
1366 #else
1367 call rf_self_excludes(i, j, sw, 1.0_dp, eFrame, d_atm, rVal, &
1368 vpair, pot(ELECTROSTATIC_POT), f, t, do_pot)
1369 #endif
1370 endif
1371 enddo
1372 endif
1373
1374 if (do_box_dipole) then
1375 #ifdef IS_MPI
1376 call accumulate_box_dipole(i, eFrame, q(:,i), pChg_local, &
1377 nChg_local, pChgPos_local, nChgPos_local, dipVec_local, &
1378 pChgCount_local, nChgCount_local)
1379 #else
1380 call accumulate_box_dipole(i, eFrame, q(:,i), pChg, nChg, &
1381 pChgPos, nChgPos, dipVec, pChgCount, nChgCount)
1382 #endif
1383 endif
1384 enddo
1385 endif
1386
1387 #ifdef IS_MPI
1388 if (do_pot) then
1389 #ifdef SINGLE_PRECISION
1390 call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_real,mpi_sum, &
1391 mpi_comm_world,mpi_err)
1392 #else
1393 call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_double_precision, &
1394 mpi_sum, mpi_comm_world,mpi_err)
1395 #endif
1396 endif
1397
1398 if (do_box_dipole) then
1399
1400 #ifdef SINGLE_PRECISION
1401 call mpi_allreduce(pChg_local, pChg, 1, mpi_real, mpi_sum, &
1402 mpi_comm_world, mpi_err)
1403 call mpi_allreduce(nChg_local, nChg, 1, mpi_real, mpi_sum, &
1404 mpi_comm_world, mpi_err)
1405 call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer, mpi_sum,&
1406 mpi_comm_world, mpi_err)
1407 call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer, mpi_sum,&
1408 mpi_comm_world, mpi_err)
1409 call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_real, mpi_sum, &
1410 mpi_comm_world, mpi_err)
1411 call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_real, mpi_sum, &
1412 mpi_comm_world, mpi_err)
1413 call mpi_allreduce(dipVec_local, dipVec, 3, mpi_real, mpi_sum, &
1414 mpi_comm_world, mpi_err)
1415 #else
1416 call mpi_allreduce(pChg_local, pChg, 1, mpi_double_precision, mpi_sum, &
1417 mpi_comm_world, mpi_err)
1418 call mpi_allreduce(nChg_local, nChg, 1, mpi_double_precision, mpi_sum, &
1419 mpi_comm_world, mpi_err)
1420 call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer,&
1421 mpi_sum, mpi_comm_world, mpi_err)
1422 call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer,&
1423 mpi_sum, mpi_comm_world, mpi_err)
1424 call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_double_precision, &
1425 mpi_sum, mpi_comm_world, mpi_err)
1426 call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_double_precision, &
1427 mpi_sum, mpi_comm_world, mpi_err)
1428 call mpi_allreduce(dipVec_local, dipVec, 3, mpi_double_precision, &
1429 mpi_sum, mpi_comm_world, mpi_err)
1430 #endif
1431
1432 endif
1433
1434 #endif
1435
1436 if (do_box_dipole) then
1437 ! first load the accumulated dipole moment (if dipoles were present)
1438 boxDipole(1) = dipVec(1)
1439 boxDipole(2) = dipVec(2)
1440 boxDipole(3) = dipVec(3)
1441
1442 ! now include the dipole moment due to charges
1443 ! use the lesser of the positive and negative charge totals
1444 if (nChg .le. pChg) then
1445 chg_value = nChg
1446 else
1447 chg_value = pChg
1448 endif
1449
1450 ! find the average positions
1451 if (pChgCount .gt. 0 .and. nChgCount .gt. 0) then
1452 pChgPos = pChgPos / pChgCount
1453 nChgPos = nChgPos / nChgCount
1454 endif
1455
1456 ! dipole is from the negative to the positive (physics notation)
1457 chgVec(1) = pChgPos(1) - nChgPos(1)
1458 chgVec(2) = pChgPos(2) - nChgPos(2)
1459 chgVec(3) = pChgPos(3) - nChgPos(3)
1460
1461 boxDipole(1) = boxDipole(1) + chgVec(1) * chg_value
1462 boxDipole(2) = boxDipole(2) + chgVec(2) * chg_value
1463 boxDipole(3) = boxDipole(3) + chgVec(3) * chg_value
1464
1465 endif
1466
1467 end subroutine do_force_loop
1468
1469 subroutine do_pair(i, j, rijsq, d, sw, do_pot, &
1470 eFrame, A, f, t, pot, particle_pot, vpair, &
1471 fpair, d_grp, r_grp, rCut, topoDist)
1472
1473 real( kind = dp ) :: vpair, sw
1474 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1475 real( kind = dp ), dimension(nLocal) :: particle_pot
1476 real( kind = dp ), dimension(3) :: fpair
1477 real( kind = dp ), dimension(nLocal) :: mfact
1478 real( kind = dp ), dimension(9,nLocal) :: eFrame
1479 real( kind = dp ), dimension(9,nLocal) :: A
1480 real( kind = dp ), dimension(3,nLocal) :: f
1481 real( kind = dp ), dimension(3,nLocal) :: t
1482
1483 logical, intent(inout) :: do_pot
1484 integer, intent(in) :: i, j
1485 real ( kind = dp ), intent(inout) :: rijsq
1486 real ( kind = dp ), intent(inout) :: r_grp
1487 real ( kind = dp ), intent(inout) :: d(3)
1488 real ( kind = dp ), intent(inout) :: d_grp(3)
1489 real ( kind = dp ), intent(inout) :: rCut
1490 integer, intent(inout) :: topoDist
1491 real ( kind = dp ) :: r, pair_pot, vdwMult, electroMult
1492 real ( kind = dp ) :: a_k, b_k, c_k, d_k, dx
1493
1494 real( kind = dp), dimension(3) :: f1, t1, t2
1495 real( kind = dp), dimension(9) :: A1, A2, eF1, eF2
1496 real( kind = dp) :: dfrhodrho_i, dfrhodrho_j
1497 real( kind = dp) :: rho_i, rho_j
1498 real( kind = dp) :: fshift_i, fshift_j
1499 real( kind = dp) :: p_vdw, p_elect, p_hb, p_met
1500 integer :: atid_i, atid_j, id1, id2, idx
1501 integer :: k
1502
1503 integer :: iHash
1504
1505 !!$ write(*,*) i, j, rijsq, d(1), d(2), d(3), sw, do_pot
1506 !!$ write(*,*) particle_pot(1), vpair, fpair(1), fpair(2), fpair(3)
1507 !!$ write(*,*) rCut
1508
1509 r = sqrt(rijsq)
1510
1511 vpair = 0.0_dp
1512 fpair(1:3) = 0.0_dp
1513
1514 p_vdw = 0.0
1515 p_elect = 0.0
1516 p_hb = 0.0
1517 p_met = 0.0
1518
1519 f1(1:3) = 0.0
1520 t1(1:3) = 0.0
1521 t2(1:3) = 0.0
1522
1523 #ifdef IS_MPI
1524 atid_i = atid_row(i)
1525 atid_j = atid_col(j)
1526
1527 do idx = 1, 9
1528 A1(idx) = A_Row(idx, i)
1529 A2(idx) = A_Col(idx, j)
1530 eF1(idx) = eFrame_Row(idx, i)
1531 eF2(idx) = eFrame_Col(idx, j)
1532 enddo
1533
1534 #else
1535 atid_i = atid(i)
1536 atid_j = atid(j)
1537 do idx = 1, 9
1538 A1(idx) = A(idx, i)
1539 A2(idx) = A(idx, j)
1540 eF1(idx) = eFrame(idx, i)
1541 eF2(idx) = eFrame(idx, j)
1542 enddo
1543
1544 #endif
1545
1546
1547 iHash = InteractionHash(atid_i, atid_j)
1548
1549 !! For the metallic potentials, we need to pass dF[rho]/drho since
1550 !! the pair calculation routines no longer are aware of parallel.
1551
1552 if ( (iand(iHash, EAM_PAIR).ne.0) .or. (iand(iHash, SC_PAIR).ne.0) ) then
1553 #ifdef IS_MPI
1554 dfrhodrho_i = dfrhodrho_row(i)
1555 dfrhodrho_j = dfrhodrho_col(j)
1556 rho_i = rho_row(i)
1557 rho_j = rho_col(j)
1558 #else
1559 dfrhodrho_i = dfrhodrho(i)
1560 dfrhodrho_j = dfrhodrho(j)
1561 rho_i = rho(i)
1562 rho_j = rho(j)
1563 #endif
1564 end if
1565
1566 vdwMult = vdwScale(topoDist)
1567 electroMult = electrostaticScale(topoDist)
1568
1569 if ( iand(iHash, LJ_PAIR).ne.0 ) then
1570 call do_lj_pair(i, j, atid_i, atid_j, d, r, rijsq, rcut, sw, vdwMult, vpair, fpair, &
1571 p_vdw, f1, do_pot)
1572 endif
1573
1574 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1575 call doElectrostaticPair(i, j, atid_i, atid_j, d, r, rijsq, rcut, sw, electroMult, &
1576 vpair, fpair, p_elect, eF1, eF2, f1, t1, t2, do_pot)
1577 endif
1578
1579 if ( iand(iHash, STICKY_PAIR).ne.0 ) then
1580 call do_sticky_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1581 p_hb, A1, A2, f1, t1, t2, do_pot)
1582 endif
1583
1584 if ( iand(iHash, STICKYPOWER_PAIR).ne.0 ) then
1585 call do_sticky_power_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1586 p_hb, A1, A2, f1, t1, t2, do_pot)
1587 endif
1588
1589 if ( iand(iHash, GAYBERNE_PAIR).ne.0 ) then
1590 call do_gb_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vdwMult, vpair, fpair, &
1591 p_vdw, A1, A2, f1, t1, t2, do_pot)
1592 endif
1593
1594 if ( iand(iHash, GAYBERNE_LJ).ne.0 ) then
1595 call do_gb_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vdwMult, vpair, fpair, &
1596 p_vdw, A1, A2, f1, t1, t2, do_pot)
1597 endif
1598
1599 if ( iand(iHash, SHAPE_PAIR).ne.0 ) then
1600 call do_shape_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1601 p_vdw, A1, A2, f1, t1, t2, do_pot)
1602 endif
1603
1604 if ( iand(iHash, SHAPE_LJ).ne.0 ) then
1605 call do_shape_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1606 p_vdw, A1, A2, f1, t1, t2, do_pot)
1607 endif
1608
1609 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1610 call do_eam_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vpair, &
1611 fpair, p_met, f1, rho_i, rho_j, dfrhodrho_i, dfrhodrho_j, fshift_i,fshift_j, do_pot)
1612 endif
1613
1614 if ( iand(iHash, SC_PAIR).ne.0 ) then
1615 call do_SC_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vpair, &
1616 fpair, p_met, f1, rho_i, rho_j, dfrhodrho_i, dfrhodrho_j, fshift_i, fshift_j, do_pot)
1617 endif
1618
1619 if ( iand(iHash, MNM_PAIR).ne.0 ) then
1620 call do_mnm_pair(i, j, atid_i, atid_j, d, r, rijsq, rcut, sw, vdwMult, vpair, fpair, &
1621 p_vdw, A1, A2, f1, t1, t2, do_pot)
1622 endif
1623
1624
1625 #ifdef IS_MPI
1626 id1 = AtomRowToGlobal(i)
1627 id2 = AtomColToGlobal(j)
1628
1629 pot_row(VDW_POT,i) = pot_row(VDW_POT,i) + 0.5*p_vdw
1630 pot_col(VDW_POT,j) = pot_col(VDW_POT,j) + 0.5*p_vdw
1631 pot_row(ELECTROSTATIC_POT,i) = pot_row(ELECTROSTATIC_POT,i) + 0.5*p_elect
1632 pot_col(ELECTROSTATIC_POT,j) = pot_col(ELECTROSTATIC_POT,j) + 0.5*p_elect
1633 pot_row(HB_POT,i) = pot_row(HB_POT,i) + 0.5*p_hb
1634 pot_col(HB_POT,j) = pot_col(HB_POT,j) + 0.5*p_hb
1635 pot_Row(METALLIC_POT,i) = pot_Row(METALLIC_POT,i) + 0.5*p_met
1636 pot_Col(METALLIC_POT,j) = pot_Col(METALLIC_POT,j) + 0.5*p_met
1637
1638 do idx = 1, 3
1639 f_Row(idx,i) = f_Row(idx,i) + f1(idx)
1640 f_Col(idx,j) = f_Col(idx,j) - f1(idx)
1641
1642 t_Row(idx,i) = t_Row(idx,i) + t1(idx)
1643 t_Col(idx,j) = t_Col(idx,j) + t2(idx)
1644 enddo
1645 ! particle_pot is the difference between the full potential
1646 ! and the full potential without the presence of a particular
1647 ! particle (atom1).
1648 !
1649 ! This reduces the density at other particle locations, so
1650 ! we need to recompute the density at atom2 assuming atom1
1651 ! didn't contribute. This then requires recomputing the
1652 ! density functional for atom2 as well.
1653 !
1654 ! Most of the particle_pot heavy lifting comes from the
1655 ! pair interaction, and will be handled by vpair. Parallel version.
1656
1657 if ( (iand(iHash, EAM_PAIR).ne.0) .or. (iand(iHash, SC_PAIR).ne.0) ) then
1658 ppot_row(i) = ppot_row(i) - frho_row(j) + fshift_j
1659 ppot_col(j) = ppot_col(j) - frho_col(i) + fshift_i
1660 end if
1661
1662 #else
1663 id1 = i
1664 id2 = j
1665
1666 pot(VDW_POT) = pot(VDW_POT) + p_vdw
1667 pot(ELECTROSTATIC_POT) = pot(ELECTROSTATIC_POT) + p_elect
1668 pot(HB_POT) = pot(HB_POT) + p_hb
1669 pot(METALLIC_POT) = pot(METALLIC_POT) + p_met
1670
1671 do idx = 1, 3
1672 f(idx,i) = f(idx,i) + f1(idx)
1673 f(idx,j) = f(idx,j) - f1(idx)
1674
1675 t(idx,i) = t(idx,i) + t1(idx)
1676 t(idx,j) = t(idx,j) + t2(idx)
1677 enddo
1678 ! particle_pot is the difference between the full potential
1679 ! and the full potential without the presence of a particular
1680 ! particle (atom1).
1681 !
1682 ! This reduces the density at other particle locations, so
1683 ! we need to recompute the density at atom2 assuming atom1
1684 ! didn't contribute. This then requires recomputing the
1685 ! density functional for atom2 as well.
1686 !
1687 ! Most of the particle_pot heavy lifting comes from the
1688 ! pair interaction, and will be handled by vpair. NonParallel version.
1689
1690 if ( (iand(iHash, EAM_PAIR).ne.0) .or. (iand(iHash, SC_PAIR).ne.0) ) then
1691 particle_pot(i) = particle_pot(i) - frho(j) + fshift_j
1692 particle_pot(j) = particle_pot(j) - frho(i) + fshift_i
1693 end if
1694
1695
1696 #endif
1697
1698 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
1699
1700 fpair(1) = fpair(1) + f1(1)
1701 fpair(2) = fpair(2) + f1(2)
1702 fpair(3) = fpair(3) + f1(3)
1703
1704 endif
1705
1706
1707 !!$
1708 !!$ particle_pot(i) = particle_pot(i) + vpair*sw
1709 !!$ particle_pot(j) = particle_pot(j) + vpair*sw
1710
1711 end subroutine do_pair
1712
1713 subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, rCut, &
1714 do_pot, do_stress, eFrame, A, f, t, pot)
1715
1716 real( kind = dp ) :: sw
1717 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1718 real( kind = dp ), dimension(9,nLocal) :: eFrame
1719 real (kind=dp), dimension(9,nLocal) :: A
1720 real (kind=dp), dimension(3,nLocal) :: f
1721 real (kind=dp), dimension(3,nLocal) :: t
1722
1723 logical, intent(inout) :: do_pot, do_stress
1724 integer, intent(in) :: i, j
1725 real ( kind = dp ), intent(inout) :: rijsq, rcijsq, rCut
1726 real ( kind = dp ) :: r, rc
1727 real ( kind = dp ), intent(inout) :: d(3), dc(3)
1728 real ( kind = dp ) :: rho_i_at_j, rho_j_at_i
1729 integer :: atid_i, atid_j, iHash
1730
1731 r = sqrt(rijsq)
1732
1733 #ifdef IS_MPI
1734 atid_i = atid_row(i)
1735 atid_j = atid_col(j)
1736 #else
1737 atid_i = atid(i)
1738 atid_j = atid(j)
1739 #endif
1740 rho_i_at_j = 0.0_dp
1741 rho_j_at_i = 0.0_dp
1742
1743 iHash = InteractionHash(atid_i, atid_j)
1744
1745 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1746 call calc_EAM_prepair_rho(i, j, atid_i, atid_j, d, r, rijsq, rho_i_at_j, rho_j_at_i)
1747 endif
1748
1749 if ( iand(iHash, SC_PAIR).ne.0 ) then
1750 call calc_SC_prepair_rho(i, j, atid_i, atid_j, d, r, rijsq, rho_i_at_j, rho_j_at_i)
1751 endif
1752
1753 if ( iand(iHash, EAM_PAIR).ne.0 .or. iand(iHash, SC_PAIR).ne.0 ) then
1754 #ifdef IS_MPI
1755 rho_col(j) = rho_col(j) + rho_i_at_j
1756 rho_row(i) = rho_row(i) + rho_j_at_i
1757 #else
1758 rho(j) = rho(j) + rho_i_at_j
1759 rho(i) = rho(i) + rho_j_at_i
1760 #endif
1761 endif
1762
1763 end subroutine do_prepair
1764
1765
1766 subroutine do_preforce(nlocal, pot, particle_pot)
1767 integer :: nlocal
1768 real( kind = dp ),dimension(LR_POT_TYPES) :: pot
1769 real( kind = dp ),dimension(nlocal) :: particle_pot
1770 integer :: sc_err = 0
1771
1772 #ifdef IS_MPI
1773 if ((FF_uses_EAM .and. SIM_uses_EAM) .or. (FF_uses_SC .and. SIM_uses_SC)) then
1774 call scatter(rho_row,rho,plan_atom_row,sc_err)
1775 if (sc_err /= 0 ) then
1776 call handleError("do_preforce()", "Error scattering rho_row into rho")
1777 endif
1778 call scatter(rho_col,rho_tmp,plan_atom_col,sc_err)
1779 if (sc_err /= 0 ) then
1780 call handleError("do_preforce()", "Error scattering rho_col into rho")
1781 endif
1782 rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal)
1783 end if
1784 #endif
1785
1786
1787
1788 if (FF_uses_EAM .and. SIM_uses_EAM) then
1789 call calc_EAM_preforce_Frho(nlocal, pot(METALLIC_POT), particle_pot)
1790 endif
1791 if (FF_uses_SC .and. SIM_uses_SC) then
1792 call calc_SC_preforce_Frho(nlocal, pot(METALLIC_POT), particle_pot)
1793 endif
1794
1795 #ifdef IS_MPI
1796 if ((FF_uses_EAM .and. SIM_uses_EAM) .or. (FF_uses_SC .and. SIM_uses_SC)) then
1797 !! communicate f(rho) and derivatives back into row and column arrays
1798 call gather(frho,frho_row,plan_atom_row, sc_err)
1799 if (sc_err /= 0) then
1800 call handleError("do_preforce()","MPI gather frho_row failure")
1801 endif
1802 call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, sc_err)
1803 if (sc_err /= 0) then
1804 call handleError("do_preforce()","MPI gather dfrhodrho_row failure")
1805 endif
1806 call gather(frho,frho_col,plan_atom_col, sc_err)
1807 if (sc_err /= 0) then
1808 call handleError("do_preforce()","MPI gather frho_col failure")
1809 endif
1810 call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, sc_err)
1811 if (sc_err /= 0) then
1812 call handleError("do_preforce()","MPI gather dfrhodrho_col failure")
1813 endif
1814 end if
1815 #endif
1816
1817 end subroutine do_preforce
1818
1819
1820 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
1821
1822 real (kind = dp), dimension(3) :: q_i
1823 real (kind = dp), dimension(3) :: q_j
1824 real ( kind = dp ), intent(out) :: r_sq
1825 real( kind = dp ) :: d(3), scaled(3)
1826 integer i
1827
1828 d(1) = q_j(1) - q_i(1)
1829 d(2) = q_j(2) - q_i(2)
1830 d(3) = q_j(3) - q_i(3)
1831
1832 ! Wrap back into periodic box if necessary
1833 if ( SIM_uses_PBC ) then
1834
1835 if( .not.boxIsOrthorhombic ) then
1836 ! calc the scaled coordinates.
1837 ! scaled = matmul(HmatInv, d)
1838
1839 scaled(1) = HmatInv(1,1)*d(1) + HmatInv(1,2)*d(2) + HmatInv(1,3)*d(3)
1840 scaled(2) = HmatInv(2,1)*d(1) + HmatInv(2,2)*d(2) + HmatInv(2,3)*d(3)
1841 scaled(3) = HmatInv(3,1)*d(1) + HmatInv(3,2)*d(2) + HmatInv(3,3)*d(3)
1842
1843 ! wrap the scaled coordinates
1844
1845 scaled(1) = scaled(1) - anint(scaled(1), kind=dp)
1846 scaled(2) = scaled(2) - anint(scaled(2), kind=dp)
1847 scaled(3) = scaled(3) - anint(scaled(3), kind=dp)
1848
1849 ! calc the wrapped real coordinates from the wrapped scaled
1850 ! coordinates
1851 ! d = matmul(Hmat,scaled)
1852 d(1)= Hmat(1,1)*scaled(1) + Hmat(1,2)*scaled(2) + Hmat(1,3)*scaled(3)
1853 d(2)= Hmat(2,1)*scaled(1) + Hmat(2,2)*scaled(2) + Hmat(2,3)*scaled(3)
1854 d(3)= Hmat(3,1)*scaled(1) + Hmat(3,2)*scaled(2) + Hmat(3,3)*scaled(3)
1855
1856 else
1857 ! calc the scaled coordinates.
1858
1859 scaled(1) = d(1) * HmatInv(1,1)
1860 scaled(2) = d(2) * HmatInv(2,2)
1861 scaled(3) = d(3) * HmatInv(3,3)
1862
1863 ! wrap the scaled coordinates
1864
1865 scaled(1) = scaled(1) - anint(scaled(1), kind=dp)
1866 scaled(2) = scaled(2) - anint(scaled(2), kind=dp)
1867 scaled(3) = scaled(3) - anint(scaled(3), kind=dp)
1868
1869 ! calc the wrapped real coordinates from the wrapped scaled
1870 ! coordinates
1871
1872 d(1) = scaled(1)*Hmat(1,1)
1873 d(2) = scaled(2)*Hmat(2,2)
1874 d(3) = scaled(3)*Hmat(3,3)
1875
1876 endif
1877
1878 endif
1879
1880 r_sq = d(1)*d(1) + d(2)*d(2) + d(3)*d(3)
1881
1882 end subroutine get_interatomic_vector
1883
1884 subroutine zero_work_arrays()
1885
1886 #ifdef IS_MPI
1887
1888 q_Row = 0.0_dp
1889 q_Col = 0.0_dp
1890
1891 q_group_Row = 0.0_dp
1892 q_group_Col = 0.0_dp
1893
1894 eFrame_Row = 0.0_dp
1895 eFrame_Col = 0.0_dp
1896
1897 A_Row = 0.0_dp
1898 A_Col = 0.0_dp
1899
1900 f_Row = 0.0_dp
1901 f_Col = 0.0_dp
1902 f_Temp = 0.0_dp
1903
1904 t_Row = 0.0_dp
1905 t_Col = 0.0_dp
1906 t_Temp = 0.0_dp
1907
1908 pot_Row = 0.0_dp
1909 pot_Col = 0.0_dp
1910 pot_Temp = 0.0_dp
1911 ppot_Temp = 0.0_dp
1912
1913 frho_row = 0.0_dp
1914 frho_col = 0.0_dp
1915 rho_row = 0.0_dp
1916 rho_col = 0.0_dp
1917 rho_tmp = 0.0_dp
1918 dfrhodrho_row = 0.0_dp
1919 dfrhodrho_col = 0.0_dp
1920
1921 #endif
1922 rho = 0.0_dp
1923 frho = 0.0_dp
1924 dfrhodrho = 0.0_dp
1925
1926 end subroutine zero_work_arrays
1927
1928 function skipThisPair(atom1, atom2) result(skip_it)
1929 integer, intent(in) :: atom1
1930 integer, intent(in), optional :: atom2
1931 logical :: skip_it
1932 integer :: unique_id_1, unique_id_2
1933 integer :: me_i,me_j
1934 integer :: i
1935
1936 skip_it = .false.
1937
1938 !! there are a number of reasons to skip a pair or a particle
1939 !! mostly we do this to exclude atoms who are involved in short
1940 !! range interactions (bonds, bends, torsions), but we also need
1941 !! to exclude some overcounted interactions that result from
1942 !! the parallel decomposition
1943
1944 #ifdef IS_MPI
1945 !! in MPI, we have to look up the unique IDs for each atom
1946 unique_id_1 = AtomRowToGlobal(atom1)
1947 unique_id_2 = AtomColToGlobal(atom2)
1948 !! this situation should only arise in MPI simulations
1949 if (unique_id_1 == unique_id_2) then
1950 skip_it = .true.
1951 return
1952 end if
1953
1954 !! this prevents us from doing the pair on multiple processors
1955 if (unique_id_1 < unique_id_2) then
1956 if (mod(unique_id_1 + unique_id_2,2) == 0) then
1957 skip_it = .true.
1958 return
1959 endif
1960 else
1961 if (mod(unique_id_1 + unique_id_2,2) == 1) then
1962 skip_it = .true.
1963 return
1964 endif
1965 endif
1966 #else
1967 !! in the normal loop, the atom numbers are unique
1968 unique_id_1 = atom1
1969 unique_id_2 = atom2
1970 #endif
1971
1972 #ifdef IS_MPI
1973 do i = 1, nSkipsForRowAtom(atom1)
1974 if (skipsForRowAtom(atom1, i) .eq. unique_id_2) then
1975 skip_it = .true.
1976 return
1977 endif
1978 end do
1979 #else
1980 do i = 1, nSkipsForLocalAtom(atom1)
1981 if (skipsForLocalAtom(atom1, i) .eq. unique_id_2) then
1982 skip_it = .true.
1983 return
1984 endif
1985 end do
1986 #endif
1987
1988 return
1989 end function skipThisPair
1990
1991 function getTopoDistance(atom1, atom2) result(topoDist)
1992 integer, intent(in) :: atom1
1993 integer, intent(in) :: atom2
1994 integer :: topoDist
1995 integer :: unique_id_2
1996 integer :: i
1997
1998 #ifdef IS_MPI
1999 unique_id_2 = AtomColToGlobal(atom2)
2000 #else
2001 unique_id_2 = atom2
2002 #endif
2003
2004 ! zero is default for unconnected (i.e. normal) pair interactions
2005
2006 topoDist = 0
2007
2008 do i = 1, nTopoPairsForAtom(atom1)
2009 if (toposForAtom(atom1, i) .eq. unique_id_2) then
2010 topoDist = topoDistance(atom1, i)
2011 return
2012 endif
2013 end do
2014
2015 return
2016 end function getTopoDistance
2017
2018 function FF_UsesDirectionalAtoms() result(doesit)
2019 logical :: doesit
2020 doesit = FF_uses_DirectionalAtoms
2021 end function FF_UsesDirectionalAtoms
2022
2023 function FF_RequiresPrepairCalc() result(doesit)
2024 logical :: doesit
2025 doesit = FF_uses_EAM .or. FF_uses_SC
2026 end function FF_RequiresPrepairCalc
2027
2028 #ifdef PROFILE
2029 function getforcetime() result(totalforcetime)
2030 real(kind=dp) :: totalforcetime
2031 totalforcetime = forcetime
2032 end function getforcetime
2033 #endif
2034
2035 !! This cleans componets of force arrays belonging only to fortran
2036
2037 subroutine add_stress_tensor(dpair, fpair, tau)
2038
2039 real( kind = dp ), dimension(3), intent(in) :: dpair, fpair
2040 real( kind = dp ), dimension(9), intent(inout) :: tau
2041
2042 ! because the d vector is the rj - ri vector, and
2043 ! because fx, fy, fz are the force on atom i, we need a
2044 ! negative sign here:
2045
2046 tau(1) = tau(1) - dpair(1) * fpair(1)
2047 tau(2) = tau(2) - dpair(1) * fpair(2)
2048 tau(3) = tau(3) - dpair(1) * fpair(3)
2049 tau(4) = tau(4) - dpair(2) * fpair(1)
2050 tau(5) = tau(5) - dpair(2) * fpair(2)
2051 tau(6) = tau(6) - dpair(2) * fpair(3)
2052 tau(7) = tau(7) - dpair(3) * fpair(1)
2053 tau(8) = tau(8) - dpair(3) * fpair(2)
2054 tau(9) = tau(9) - dpair(3) * fpair(3)
2055
2056 end subroutine add_stress_tensor
2057
2058 end module doForces