ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/heatflux/src/UseTheForce/doForces.F90
Revision: 1685
Committed: Thu Mar 1 21:22:42 2012 UTC (13 years, 5 months ago) by chuckv
File size: 70567 byte(s)
Log Message:
Bug squashes. Pre-existing bug where ppot_row and ppot_col were never zeroed. In non-metallic potentials these arrays should be zero.

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

Properties

Name Value
svn:keywords Author Id Revision Date