ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/doForces.F90
Revision: 1689
Committed: Wed Mar 14 17:57:08 2012 UTC (13 years, 1 month ago) by gezelter
File size: 68872 byte(s)
Log Message:
Bug fixes for GB.  Now using strength parameter mixing ideas from Wu
et al. [J. Chem. Phys. 135, 155104 (2011)].  This helps get the
dissimilar particle mixing behavior to be the same whichever order the
two particles come in.  This does require that the force field file to
specify explicitly the values for epsilon in the cross (X), side-by-side (S), 
and end-to-end (E) configurations.


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

Properties

Name Value
svn:keywords Author Id Revision Date