ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/doForces.F90
Revision: 1615
Committed: Fri Aug 26 17:55:44 2011 UTC (13 years, 8 months ago) by gezelter
File size: 68628 byte(s)
Log Message:
Added Momentum correlation function, imported changes from Vector from
development branch, updated comments in some integrators

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 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
1101 eFrame, A, f, t, pot, particle_pot, vpair, &
1102 fpair, d_grp, rgrp, rCut, topoDist)
1103 #endif
1104 vij = vij + vpair
1105 fij(1) = fij(1) + fpair(1)
1106 fij(2) = fij(2) + fpair(2)
1107 fij(3) = fij(3) + fpair(3)
1108 call add_stress_tensor(d_atm, fpair, tau)
1109 endif
1110 enddo inner
1111 enddo
1112
1113 if (loop .eq. PAIR_LOOP) then
1114 if (in_switching_region) then
1115 swderiv = vij*dswdr/rgrp
1116 fg = swderiv*d_grp
1117 fij(1) = fij(1) + fg(1)
1118 fij(2) = fij(2) + fg(2)
1119 fij(3) = fij(3) + fg(3)
1120
1121 if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
1122 call add_stress_tensor(d_atm, fg, tau)
1123 endif
1124
1125 do ia=groupStartRow(i), groupStartRow(i+1)-1
1126 atom1=groupListRow(ia)
1127 mf = mfactRow(atom1)
1128 ! fg is the force on atom ia due to cutoff group's
1129 ! presence in switching region
1130 fg = swderiv*d_grp*mf
1131 #ifdef IS_MPI
1132 f_Row(1,atom1) = f_Row(1,atom1) + fg(1)
1133 f_Row(2,atom1) = f_Row(2,atom1) + fg(2)
1134 f_Row(3,atom1) = f_Row(3,atom1) + fg(3)
1135 #else
1136 f(1,atom1) = f(1,atom1) + fg(1)
1137 f(2,atom1) = f(2,atom1) + fg(2)
1138 f(3,atom1) = f(3,atom1) + fg(3)
1139 #endif
1140 if (n_in_i .gt. 1) then
1141 if (SIM_uses_AtomicVirial) then
1142 ! find the distance between the atom
1143 ! and the center of the cutoff group:
1144 #ifdef IS_MPI
1145 call get_interatomic_vector(q_Row(:,atom1), &
1146 q_group_Row(:,i), dag, rag)
1147 #else
1148 call get_interatomic_vector(q(:,atom1), &
1149 q_group(:,i), dag, rag)
1150 #endif
1151 call add_stress_tensor(dag,fg,tau)
1152 endif
1153 endif
1154 enddo
1155
1156 do jb=groupStartCol(j), groupStartCol(j+1)-1
1157 atom2=groupListCol(jb)
1158 mf = mfactCol(atom2)
1159 ! fg is the force on atom jb due to cutoff group's
1160 ! presence in switching region
1161 fg = -swderiv*d_grp*mf
1162 #ifdef IS_MPI
1163 f_Col(1,atom2) = f_Col(1,atom2) + fg(1)
1164 f_Col(2,atom2) = f_Col(2,atom2) + fg(2)
1165 f_Col(3,atom2) = f_Col(3,atom2) + fg(3)
1166 #else
1167 f(1,atom2) = f(1,atom2) + fg(1)
1168 f(2,atom2) = f(2,atom2) + fg(2)
1169 f(3,atom2) = f(3,atom2) + fg(3)
1170 #endif
1171 if (n_in_j .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_Col(:,atom2), &
1177 q_group_Col(:,j), dag, rag)
1178 #else
1179 call get_interatomic_vector(q(:,atom2), &
1180 q_group(:,j), dag, rag)
1181 #endif
1182 call add_stress_tensor(dag,fg,tau)
1183 endif
1184 endif
1185 enddo
1186 endif
1187 !if (.not.SIM_uses_AtomicVirial) then
1188 ! call add_stress_tensor(d_grp, fij, tau)
1189 !endif
1190 endif
1191 endif
1192 endif
1193 enddo
1194
1195 enddo outer
1196
1197 if (update_nlist) then
1198 #ifdef IS_MPI
1199 point(nGroupsInRow + 1) = nlist + 1
1200 #else
1201 point(nGroups) = nlist + 1
1202 #endif
1203 if (loop .eq. PREPAIR_LOOP) then
1204 ! we just did the neighbor list update on the first
1205 ! pass, so we don't need to do it
1206 ! again on the second pass
1207 update_nlist = .false.
1208 endif
1209 endif
1210
1211 if (loop .eq. PREPAIR_LOOP) then
1212 #ifdef IS_MPI
1213 call do_preforce(nlocal, pot_local, particle_pot)
1214 #else
1215 call do_preforce(nlocal, pot, particle_pot)
1216 #endif
1217 endif
1218
1219 enddo
1220
1221 !! Do timing
1222 #ifdef PROFILE
1223 call cpu_time(forceTimeFinal)
1224 forceTime = forceTime + forceTimeFinal - forceTimeInitial
1225 #endif
1226
1227 #ifdef IS_MPI
1228 !!distribute forces
1229
1230 f_temp = 0.0_dp
1231 call scatter(f_Row,f_temp,plan_atom_row_3d)
1232 do i = 1,nlocal
1233 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1234 end do
1235
1236 f_temp = 0.0_dp
1237 call scatter(f_Col,f_temp,plan_atom_col_3d)
1238 do i = 1,nlocal
1239 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1240 end do
1241
1242 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
1243 t_temp = 0.0_dp
1244 call scatter(t_Row,t_temp,plan_atom_row_3d)
1245 do i = 1,nlocal
1246 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1247 end do
1248 t_temp = 0.0_dp
1249 call scatter(t_Col,t_temp,plan_atom_col_3d)
1250
1251 do i = 1,nlocal
1252 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1253 end do
1254 endif
1255
1256 ! scatter/gather pot_row into the members of my column
1257 do i = 1,LR_POT_TYPES
1258 call scatter(pot_Row(i,:), pot_Temp(i,:), plan_atom_row)
1259 end do
1260 ! scatter/gather pot_local into all other procs
1261 ! add resultant to get total pot
1262 do i = 1, nlocal
1263 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES) &
1264 + pot_Temp(1:LR_POT_TYPES,i)
1265 enddo
1266
1267 ! factor of two is because the total potential terms are divided by 2 in parallel
1268 ! due to row/ column scatter
1269 do i = 1,LR_POT_TYPES
1270 particle_pot(1:nlocal) = particle_pot(1:nlocal) + 2.0 * pot_Temp(i,1:nlocal)
1271 enddo
1272
1273
1274 pot_Temp = 0.0_DP
1275
1276 do i = 1,LR_POT_TYPES
1277 call scatter(pot_Col(i,:), pot_Temp(i,:), plan_atom_col)
1278 end do
1279
1280 do i = 1, nlocal
1281 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES)&
1282 + pot_Temp(1:LR_POT_TYPES,i)
1283 enddo
1284
1285 ! factor of two is because the total potential terms are divided by 2 in parallel
1286 ! due to row/ column scatter
1287 do i = 1,LR_POT_TYPES
1288 particle_pot(1:nlocal) = particle_pot(1:nlocal) + 2.0 * pot_Temp(i,1:nlocal)
1289 enddo
1290
1291 ppot_Temp = 0.0_DP
1292
1293 call scatter(ppot_Row(:), ppot_Temp(:), plan_atom_row)
1294 do i = 1, nlocal
1295 particle_pot(i) = particle_pot(i) + ppot_Temp(i)
1296 enddo
1297
1298 ppot_Temp = 0.0_DP
1299
1300 call scatter(ppot_Col(:), ppot_Temp(:), plan_atom_col)
1301 do i = 1, nlocal
1302 particle_pot(i) = particle_pot(i) + ppot_Temp(i)
1303 enddo
1304
1305 #endif
1306
1307 if (SIM_requires_postpair_calc) then
1308 do i = 1, nlocal
1309
1310 ! we loop only over the local atoms, so we don't need row and column
1311 ! lookups for the types
1312
1313 me_i = atid(i)
1314
1315 ! is the atom electrostatic? See if it would have an
1316 ! electrostatic interaction with itself
1317 iHash = InteractionHash(me_i,me_i)
1318
1319 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1320
1321 ! loop over the excludes to accumulate charge in the
1322 ! cutoff sphere that we've left out of the normal pair loop
1323 skch = 0.0_dp
1324
1325 do i1 = 1, nSkipsForLocalAtom(i)
1326 j = skipsForLocalAtom(i, i1)
1327 me_j = atid(j)
1328 jHash = InteractionHash(me_i,me_j)
1329 if ( iand(jHash, ELECTROSTATIC_PAIR).ne.0 ) then
1330 skch = skch + getCharge(me_j)
1331 endif
1332 enddo
1333
1334 #ifdef IS_MPI
1335 call self_self(i, eFrame, skch, pot_local(ELECTROSTATIC_POT), t)
1336 #else
1337 call self_self(i, eFrame, skch, pot(ELECTROSTATIC_POT), t)
1338 #endif
1339 endif
1340
1341
1342 if (electrostaticSummationMethod.eq.REACTION_FIELD) then
1343
1344 ! loop over the excludes to accumulate RF stuff we've
1345 ! left out of the normal pair loop
1346
1347 do i1 = 1, nSkipsForLocalAtom(i)
1348 j = skipsForLocalAtom(i, i1)
1349
1350 ! prevent overcounting of the skips
1351 if (i.lt.j) then
1352 call get_interatomic_vector(q(:,i), q(:,j), d_atm, ratmsq)
1353 rVal = sqrt(ratmsq)
1354 call get_switch(ratmsq, sw, dswdr, rVal,in_switching_region)
1355 #ifdef IS_MPI
1356 call rf_self_excludes(i, j, sw, 1.0_dp, eFrame, d_atm, rVal, &
1357 vpair, pot_local(ELECTROSTATIC_POT), f, t)
1358 #else
1359 call rf_self_excludes(i, j, sw, 1.0_dp, eFrame, d_atm, rVal, &
1360 vpair, pot(ELECTROSTATIC_POT), f, t)
1361 #endif
1362 endif
1363 enddo
1364 endif
1365
1366 if (do_box_dipole) then
1367 #ifdef IS_MPI
1368 call accumulate_box_dipole(i, eFrame, q(:,i), pChg_local, &
1369 nChg_local, pChgPos_local, nChgPos_local, dipVec_local, &
1370 pChgCount_local, nChgCount_local)
1371 #else
1372 call accumulate_box_dipole(i, eFrame, q(:,i), pChg, nChg, &
1373 pChgPos, nChgPos, dipVec, pChgCount, nChgCount)
1374 #endif
1375 endif
1376 enddo
1377 endif
1378
1379 #ifdef IS_MPI
1380 #ifdef SINGLE_PRECISION
1381 call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_real,mpi_sum, &
1382 mpi_comm_world,mpi_err)
1383 #else
1384 call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_double_precision, &
1385 mpi_sum, mpi_comm_world,mpi_err)
1386 #endif
1387
1388 if (do_box_dipole) then
1389
1390 #ifdef SINGLE_PRECISION
1391 call mpi_allreduce(pChg_local, pChg, 1, mpi_real, mpi_sum, &
1392 mpi_comm_world, mpi_err)
1393 call mpi_allreduce(nChg_local, nChg, 1, mpi_real, mpi_sum, &
1394 mpi_comm_world, mpi_err)
1395 call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer, mpi_sum,&
1396 mpi_comm_world, mpi_err)
1397 call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer, mpi_sum,&
1398 mpi_comm_world, mpi_err)
1399 call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_real, mpi_sum, &
1400 mpi_comm_world, mpi_err)
1401 call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_real, mpi_sum, &
1402 mpi_comm_world, mpi_err)
1403 call mpi_allreduce(dipVec_local, dipVec, 3, mpi_real, mpi_sum, &
1404 mpi_comm_world, mpi_err)
1405 #else
1406 call mpi_allreduce(pChg_local, pChg, 1, mpi_double_precision, mpi_sum, &
1407 mpi_comm_world, mpi_err)
1408 call mpi_allreduce(nChg_local, nChg, 1, mpi_double_precision, mpi_sum, &
1409 mpi_comm_world, mpi_err)
1410 call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer,&
1411 mpi_sum, mpi_comm_world, mpi_err)
1412 call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer,&
1413 mpi_sum, mpi_comm_world, mpi_err)
1414 call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_double_precision, &
1415 mpi_sum, mpi_comm_world, mpi_err)
1416 call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_double_precision, &
1417 mpi_sum, mpi_comm_world, mpi_err)
1418 call mpi_allreduce(dipVec_local, dipVec, 3, mpi_double_precision, &
1419 mpi_sum, mpi_comm_world, mpi_err)
1420 #endif
1421
1422 endif
1423
1424 #endif
1425
1426 if (do_box_dipole) then
1427 ! first load the accumulated dipole moment (if dipoles were present)
1428 boxDipole(1) = dipVec(1)
1429 boxDipole(2) = dipVec(2)
1430 boxDipole(3) = dipVec(3)
1431
1432 ! now include the dipole moment due to charges
1433 ! use the lesser of the positive and negative charge totals
1434 if (nChg .le. pChg) then
1435 chg_value = nChg
1436 else
1437 chg_value = pChg
1438 endif
1439
1440 ! find the average positions
1441 if (pChgCount .gt. 0 .and. nChgCount .gt. 0) then
1442 pChgPos = pChgPos / pChgCount
1443 nChgPos = nChgPos / nChgCount
1444 endif
1445
1446 ! dipole is from the negative to the positive (physics notation)
1447 chgVec(1) = pChgPos(1) - nChgPos(1)
1448 chgVec(2) = pChgPos(2) - nChgPos(2)
1449 chgVec(3) = pChgPos(3) - nChgPos(3)
1450
1451 boxDipole(1) = boxDipole(1) + chgVec(1) * chg_value
1452 boxDipole(2) = boxDipole(2) + chgVec(2) * chg_value
1453 boxDipole(3) = boxDipole(3) + chgVec(3) * chg_value
1454
1455 endif
1456
1457 end subroutine do_force_loop
1458
1459 subroutine do_pair(i, j, rijsq, d, sw, &
1460 eFrame, A, f, t, pot, particle_pot, vpair, &
1461 fpair, d_grp, r_grp, rCut, topoDist)
1462
1463 real( kind = dp ) :: vpair, sw
1464 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1465 real( kind = dp ), dimension(nLocal) :: particle_pot
1466 real( kind = dp ), dimension(3) :: fpair
1467 real( kind = dp ), dimension(nLocal) :: mfact
1468 real( kind = dp ), dimension(9,nLocal) :: eFrame
1469 real( kind = dp ), dimension(9,nLocal) :: A
1470 real( kind = dp ), dimension(3,nLocal) :: f
1471 real( kind = dp ), dimension(3,nLocal) :: t
1472
1473 integer, intent(in) :: i, j
1474 real ( kind = dp ), intent(inout) :: rijsq
1475 real ( kind = dp ), intent(inout) :: r_grp
1476 real ( kind = dp ), intent(inout) :: d(3)
1477 real ( kind = dp ), intent(inout) :: d_grp(3)
1478 real ( kind = dp ), intent(inout) :: rCut
1479 integer, intent(inout) :: topoDist
1480 real ( kind = dp ) :: r, pair_pot, vdwMult, electroMult
1481 real ( kind = dp ) :: a_k, b_k, c_k, d_k, dx
1482
1483 real( kind = dp), dimension(3) :: f1, t1, t2
1484 real( kind = dp), dimension(9) :: A1, A2, eF1, eF2
1485 real( kind = dp) :: dfrhodrho_i, dfrhodrho_j
1486 real( kind = dp) :: rho_i, rho_j
1487 real( kind = dp) :: fshift_i, fshift_j
1488 real( kind = dp) :: p_vdw, p_elect, p_hb, p_met
1489 integer :: atid_i, atid_j, id1, id2, idx
1490 integer :: k
1491
1492 integer :: iHash
1493
1494 r = sqrt(rijsq)
1495
1496 vpair = 0.0_dp
1497 fpair(1:3) = 0.0_dp
1498
1499 p_vdw = 0.0
1500 p_elect = 0.0
1501 p_hb = 0.0
1502 p_met = 0.0
1503
1504 f1(1:3) = 0.0
1505
1506 #ifdef IS_MPI
1507 atid_i = atid_row(i)
1508 atid_j = atid_col(j)
1509 #else
1510 atid_i = atid(i)
1511 atid_j = atid(j)
1512 #endif
1513
1514 iHash = InteractionHash(atid_i, atid_j)
1515
1516 vdwMult = vdwScale(topoDist)
1517 electroMult = electrostaticScale(topoDist)
1518
1519 if ( iand(iHash, LJ_PAIR).ne.0 ) then
1520 call do_lj_pair(atid_i, atid_j, d, r, rijsq, rcut, sw, vdwMult, vpair, fpair, &
1521 p_vdw, f1)
1522 endif
1523
1524 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1525 #ifdef IS_MPI
1526 call doElectrostaticPair(atid_i, atid_j, d, r, rijsq, rcut, sw, electroMult, &
1527 vpair, fpair, p_elect, eFrame_Row(:,i), eFrame_Col(:,j), &
1528 f1, t_Row(:,i), t_Col(:,j))
1529 #else
1530 call doElectrostaticPair(atid_i, atid_j, d, r, rijsq, rcut, sw, electroMult, &
1531 vpair, fpair, p_elect, eFrame(:,i), eFrame(:,j), f1, t(:,i), t(:,j))
1532 #endif
1533 endif
1534
1535 if ( iand(iHash, STICKY_PAIR).ne.0 ) then
1536 #ifdef IS_MPI
1537 call do_sticky_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1538 p_hb, A_Row(:,i), A_Col(:,j), f1, t_Row(:,i), t_Col(:,j))
1539 #else
1540 call do_sticky_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1541 p_hb, A(:,i), A(:,j), f1, t(:,i), t(:,j))
1542 #endif
1543 endif
1544
1545 if ( iand(iHash, STICKYPOWER_PAIR).ne.0 ) then
1546 #ifdef IS_MPI
1547 call do_sticky_power_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1548 p_hb, A_Row(:,i), A_Col(:,j), f1, t_Row(:,i), t_Col(:,j))
1549 #else
1550 call do_sticky_power_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1551 p_hb, A(:,i), A(:,j), f1, t(:,i), t(:,j))
1552 #endif
1553 endif
1554
1555 if ( iand(iHash, GAYBERNE_PAIR).ne.0 ) then
1556 #ifdef IS_MPI
1557 call do_gb_pair(atid_i, atid_j, d, r, rijsq, sw, vdwMult, vpair, fpair, &
1558 p_vdw, A_Row(:,i), A_Col(:,j), f1, t_Row(:,i), t_Col(:,j))
1559 #else
1560 call do_gb_pair(atid_i, atid_j, d, r, rijsq, sw, vdwMult, vpair, fpair, &
1561 p_vdw, A(:,i), A(:,j), f1, t(:,i), t(:,j))
1562 #endif
1563 endif
1564
1565 if ( iand(iHash, GAYBERNE_LJ).ne.0 ) then
1566 #ifdef IS_MPI
1567 call do_gb_pair(atid_i, atid_j, d, r, rijsq, sw, vdwMult, vpair, fpair, &
1568 p_vdw, A_Row(:,i), A_Col(:,j), f1, t_Row(:,i), t_Col(:,j))
1569 #else
1570 call do_gb_pair(atid_i, atid_j, d, r, rijsq, sw, vdwMult, vpair, fpair, &
1571 p_vdw, A(:,i), A(:,j), f1, t(:,i), t(:,j))
1572 #endif
1573 endif
1574
1575 if ( iand(iHash, SHAPE_PAIR).ne.0 ) then
1576 #ifdef IS_MPI
1577 call do_shape_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1578 p_vdw, A_Row(:,i), A_Col(:,j), f1, t_Row(:,i), t_Col(:,j))
1579 #else
1580 call do_shape_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1581 p_vdw, A(:,i), A(:,j), f1, t(:,i), t(:,j))
1582 #endif
1583 endif
1584
1585 if ( iand(iHash, SHAPE_LJ).ne.0 ) then
1586 #ifdef IS_MPI
1587 call do_shape_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1588 p_vdw, A_Row(:,i), A_Col(:,j), f1, t_Row(:,i), t_Col(:,j))
1589 #else
1590 call do_shape_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1591 p_vdw, A(:,i), A(:,j), f1, t(:,i), t(:,j))
1592 #endif
1593 endif
1594
1595 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1596 #ifdef IS_MPI
1597 call do_eam_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, &
1598 fpair, p_met, f1, rho_row(i), rho_col(j), dfrhodrho_row(i), dfrhodrho_col(j), &
1599 fshift_i, fshift_j)
1600 #else
1601 call do_eam_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, &
1602 fpair, p_met, f1, rho(i), rho(j), dfrhodrho(i), dfrhodrho(j), fshift_i, fshift_j)
1603 #endif
1604 endif
1605
1606 if ( iand(iHash, SC_PAIR).ne.0 ) then
1607 #ifdef IS_MPI
1608 call do_SC_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, &
1609 fpair, p_met, f1, rho_row(i), rho_col(j), dfrhodrho_row(i), dfrhodrho_col(j), &
1610 fshift_i, fshift_j)
1611 #else
1612 call do_SC_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, &
1613 fpair, p_met, f1, rho(i), rho(j), dfrhodrho(i), dfrhodrho(j), fshift_i, fshift_j)
1614 #endif
1615 endif
1616
1617 if ( iand(iHash, MNM_PAIR).ne.0 ) then
1618 #ifdef IS_MPI
1619 call do_mnm_pair(atid_i, atid_j, d, r, rijsq, rcut, sw, vdwMult, vpair, fpair, &
1620 p_vdw, A_Row(:,i), A_Col(:,j), f1, t_Row(:,i), t_Col(:,j))
1621 #else
1622 call do_mnm_pair(atid_i, atid_j, d, r, rijsq, rcut, sw, vdwMult, vpair, fpair, &
1623 p_vdw, A(:,i), A(:,j), f1, t(:,i), t(:,j))
1624 #endif
1625 endif
1626
1627
1628 #ifdef IS_MPI
1629 id1 = AtomRowToGlobal(i)
1630 id2 = AtomColToGlobal(j)
1631
1632 pot_row(VDW_POT,i) = pot_row(VDW_POT,i) + 0.5*p_vdw
1633 pot_col(VDW_POT,j) = pot_col(VDW_POT,j) + 0.5*p_vdw
1634 pot_row(ELECTROSTATIC_POT,i) = pot_row(ELECTROSTATIC_POT,i) + 0.5*p_elect
1635 pot_col(ELECTROSTATIC_POT,j) = pot_col(ELECTROSTATIC_POT,j) + 0.5*p_elect
1636 pot_row(HB_POT,i) = pot_row(HB_POT,i) + 0.5*p_hb
1637 pot_col(HB_POT,j) = pot_col(HB_POT,j) + 0.5*p_hb
1638 pot_Row(METALLIC_POT,i) = pot_Row(METALLIC_POT,i) + 0.5*p_met
1639 pot_Col(METALLIC_POT,j) = pot_Col(METALLIC_POT,j) + 0.5*p_met
1640
1641 do idx = 1, 3
1642 f_Row(idx,i) = f_Row(idx,i) + f1(idx)
1643 f_Col(idx,j) = f_Col(idx,j) - f1(idx)
1644 enddo
1645 ! particle_pot is the difference between the full potential
1646 ! and the full potential without the presence of a particular
1647 ! particle (atom1).
1648 !
1649 ! This reduces the density at other particle locations, so
1650 ! we need to recompute the density at atom2 assuming atom1
1651 ! didn't contribute. This then requires recomputing the
1652 ! density functional for atom2 as well.
1653 !
1654 ! Most of the particle_pot heavy lifting comes from the
1655 ! pair interaction, and will be handled by vpair. Parallel version.
1656
1657 if ( (iand(iHash, EAM_PAIR).ne.0) .or. (iand(iHash, SC_PAIR).ne.0) ) then
1658 ppot_row(i) = ppot_row(i) - frho_row(j) + fshift_j
1659 ppot_col(j) = ppot_col(j) - frho_col(i) + fshift_i
1660 end if
1661
1662 #else
1663 id1 = i
1664 id2 = j
1665
1666 pot(VDW_POT) = pot(VDW_POT) + p_vdw
1667 pot(ELECTROSTATIC_POT) = pot(ELECTROSTATIC_POT) + p_elect
1668 pot(HB_POT) = pot(HB_POT) + p_hb
1669 pot(METALLIC_POT) = pot(METALLIC_POT) + p_met
1670
1671 ! only done for single processor. In Parallel, the particle_pot
1672 ! is constructed from the row and column potentials.
1673
1674 particle_pot(i) = particle_pot(i) + vpair*sw
1675 particle_pot(j) = particle_pot(j) + vpair*sw
1676
1677 do idx = 1, 3
1678 f(idx,i) = f(idx,i) + f1(idx)
1679 f(idx,j) = f(idx,j) - f1(idx)
1680 enddo
1681 ! particle_pot is the difference between the full potential
1682 ! and the full potential without the presence of a particular
1683 ! particle (atom1).
1684 !
1685 ! This reduces the density at other particle locations, so
1686 ! we need to recompute the density at atom2 assuming atom1
1687 ! didn't contribute. This then requires recomputing the
1688 ! density functional for atom2 as well.
1689 !
1690 ! Most of the particle_pot heavy lifting comes from the
1691 ! pair interaction, and will be handled by vpair. NonParallel version.
1692
1693 if ( (iand(iHash, EAM_PAIR).ne.0) .or. (iand(iHash, SC_PAIR).ne.0) ) then
1694 particle_pot(i) = particle_pot(i) - frho(j) + fshift_j
1695 particle_pot(j) = particle_pot(j) - frho(i) + fshift_i
1696 end if
1697
1698
1699 #endif
1700
1701 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
1702
1703 fpair(1) = fpair(1) + f1(1)
1704 fpair(2) = fpair(2) + f1(2)
1705 fpair(3) = fpair(3) + f1(3)
1706
1707 endif
1708
1709
1710 end subroutine do_pair
1711
1712 subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, rCut, &
1713 eFrame, A, f, t, pot)
1714
1715 real( kind = dp ) :: sw
1716 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1717 real( kind = dp ), dimension(9,nLocal) :: eFrame
1718 real (kind=dp), dimension(9,nLocal) :: A
1719 real (kind=dp), dimension(3,nLocal) :: f
1720 real (kind=dp), dimension(3,nLocal) :: t
1721
1722 integer, intent(in) :: i, j
1723 real ( kind = dp ), intent(inout) :: rijsq, rcijsq, rCut
1724 real ( kind = dp ) :: r, rc
1725 real ( kind = dp ), intent(inout) :: d(3), dc(3)
1726 real ( kind = dp ) :: rho_i_at_j, rho_j_at_i
1727 integer :: atid_i, atid_j, iHash
1728
1729 r = sqrt(rijsq)
1730
1731 #ifdef IS_MPI
1732 atid_i = atid_row(i)
1733 atid_j = atid_col(j)
1734 #else
1735 atid_i = atid(i)
1736 atid_j = atid(j)
1737 #endif
1738 rho_i_at_j = 0.0_dp
1739 rho_j_at_i = 0.0_dp
1740
1741 iHash = InteractionHash(atid_i, atid_j)
1742
1743 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1744 call calc_EAM_prepair_rho(atid_i, atid_j, d, r, rijsq, rho_i_at_j, rho_j_at_i)
1745 endif
1746
1747 if ( iand(iHash, SC_PAIR).ne.0 ) then
1748 call calc_SC_prepair_rho(atid_i, atid_j, d, r, rijsq, rho_i_at_j, rho_j_at_i)
1749 endif
1750
1751 if ( iand(iHash, EAM_PAIR).ne.0 .or. iand(iHash, SC_PAIR).ne.0 ) then
1752 #ifdef IS_MPI
1753 rho_col(j) = rho_col(j) + rho_i_at_j
1754 rho_row(i) = rho_row(i) + rho_j_at_i
1755 #else
1756 rho(j) = rho(j) + rho_i_at_j
1757 rho(i) = rho(i) + rho_j_at_i
1758 #endif
1759 endif
1760
1761 end subroutine do_prepair
1762
1763
1764 subroutine do_preforce(nlocal, pot, particle_pot)
1765 integer :: nlocal
1766 real( kind = dp ),dimension(LR_POT_TYPES) :: pot
1767 real( kind = dp ),dimension(nlocal) :: particle_pot
1768 integer :: sc_err = 0
1769
1770 #ifdef IS_MPI
1771 if ((FF_uses_EAM .and. SIM_uses_EAM) .or. (FF_uses_SC .and. SIM_uses_SC)) then
1772 call scatter(rho_row,rho,plan_atom_row,sc_err)
1773 if (sc_err /= 0 ) then
1774 call handleError("do_preforce()", "Error scattering rho_row into rho")
1775 endif
1776 call scatter(rho_col,rho_tmp,plan_atom_col,sc_err)
1777 if (sc_err /= 0 ) then
1778 call handleError("do_preforce()", "Error scattering rho_col into rho")
1779 endif
1780 rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal)
1781 end if
1782 #endif
1783
1784
1785
1786 if (FF_uses_EAM .and. SIM_uses_EAM) then
1787 call calc_EAM_preforce_Frho(nlocal, pot(METALLIC_POT), particle_pot)
1788 endif
1789 if (FF_uses_SC .and. SIM_uses_SC) then
1790 call calc_SC_preforce_Frho(nlocal, pot(METALLIC_POT), particle_pot)
1791 endif
1792
1793 #ifdef IS_MPI
1794 if ((FF_uses_EAM .and. SIM_uses_EAM) .or. (FF_uses_SC .and. SIM_uses_SC)) then
1795 !! communicate f(rho) and derivatives back into row and column arrays
1796 call gather(frho,frho_row,plan_atom_row, sc_err)
1797 if (sc_err /= 0) then
1798 call handleError("do_preforce()","MPI gather frho_row failure")
1799 endif
1800 call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, sc_err)
1801 if (sc_err /= 0) then
1802 call handleError("do_preforce()","MPI gather dfrhodrho_row failure")
1803 endif
1804 call gather(frho,frho_col,plan_atom_col, sc_err)
1805 if (sc_err /= 0) then
1806 call handleError("do_preforce()","MPI gather frho_col failure")
1807 endif
1808 call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, sc_err)
1809 if (sc_err /= 0) then
1810 call handleError("do_preforce()","MPI gather dfrhodrho_col failure")
1811 endif
1812 end if
1813 #endif
1814
1815 end subroutine do_preforce
1816
1817
1818 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
1819
1820 real (kind = dp), dimension(3) :: q_i
1821 real (kind = dp), dimension(3) :: q_j
1822 real ( kind = dp ), intent(out) :: r_sq
1823 real( kind = dp ) :: d(3), scaled(3)
1824 real(kind=dp)::t
1825 integer i
1826
1827 d(1) = q_j(1) - q_i(1)
1828 d(2) = q_j(2) - q_i(2)
1829 d(3) = q_j(3) - q_i(3)
1830
1831 ! Wrap back into periodic box if necessary
1832 if ( SIM_uses_PBC ) then
1833
1834 if( .not.boxIsOrthorhombic ) then
1835 ! calc the scaled coordinates.
1836 ! unwrap the matmul and do things explicitly
1837 ! scaled = matmul(HmatInv, d)
1838
1839 scaled(1) = HmatInv(1,1)*d(1) + HmatInv(1,2)*d(2) + HmatInv(1,3)*d(3)
1840 scaled(2) = HmatInv(2,1)*d(1) + HmatInv(2,2)*d(2) + HmatInv(2,3)*d(3)
1841 scaled(3) = HmatInv(3,1)*d(1) + HmatInv(3,2)*d(2) + HmatInv(3,3)*d(3)
1842
1843 ! wrap the scaled coordinates (but don't use anint for speed)
1844
1845 t = scaled(1)
1846 if (t .gt. 0.0) then
1847 scaled(1) = t - floor(t + 0.5)
1848 else
1849 scaled(1) = t - ceiling(t - 0.5)
1850 endif
1851
1852 t = scaled(2)
1853 if (t .gt. 0.0) then
1854 scaled(2) = t - floor(t + 0.5)
1855 else
1856 scaled(2) = t - ceiling(t - 0.5)
1857 endif
1858
1859 t = scaled(3)
1860 if (t .gt. 0.0) then
1861 scaled(3) = t - floor(t + 0.5)
1862 else
1863 scaled(3) = t - ceiling(t - 0.5)
1864 endif
1865
1866 ! calc the wrapped real coordinates from the wrapped scaled
1867 ! coordinates
1868 ! d = matmul(Hmat,scaled)
1869 d(1)= Hmat(1,1)*scaled(1) + Hmat(1,2)*scaled(2) + Hmat(1,3)*scaled(3)
1870 d(2)= Hmat(2,1)*scaled(1) + Hmat(2,2)*scaled(2) + Hmat(2,3)*scaled(3)
1871 d(3)= Hmat(3,1)*scaled(1) + Hmat(3,2)*scaled(2) + Hmat(3,3)*scaled(3)
1872
1873 else
1874 ! calc the scaled coordinates
1875 scaled(1) = d(1) * HmatInv(1,1)
1876 scaled(2) = d(2) * HmatInv(2,2)
1877 scaled(3) = d(3) * HmatInv(3,3)
1878
1879 ! wrap the scaled coordinates
1880
1881 t = scaled(1)
1882 if (t .gt. 0.0) then
1883 scaled(1) = t - floor(t + 0.5)
1884 else
1885 scaled(1) = t - ceiling(t - 0.5)
1886 endif
1887
1888 t = scaled(2)
1889 if (t .gt. 0.0) then
1890 scaled(2) = t - floor(t + 0.5)
1891 else
1892 scaled(2) = t - ceiling(t - 0.5)
1893 endif
1894
1895 t = scaled(3)
1896 if (t .gt. 0.0) then
1897 scaled(3) = t - floor(t + 0.5)
1898 else
1899 scaled(3) = t - ceiling(t - 0.5)
1900 endif
1901
1902 ! calc the wrapped real coordinates from the wrapped scaled
1903 ! coordinates
1904
1905 d(1) = scaled(1)*Hmat(1,1)
1906 d(2) = scaled(2)*Hmat(2,2)
1907 d(3) = scaled(3)*Hmat(3,3)
1908
1909 endif
1910
1911 endif
1912
1913 r_sq = d(1)*d(1) + d(2)*d(2) + d(3)*d(3)
1914
1915 end subroutine get_interatomic_vector
1916
1917 subroutine zero_work_arrays()
1918
1919 #ifdef IS_MPI
1920
1921 q_Row = 0.0_dp
1922 q_Col = 0.0_dp
1923
1924 q_group_Row = 0.0_dp
1925 q_group_Col = 0.0_dp
1926
1927 eFrame_Row = 0.0_dp
1928 eFrame_Col = 0.0_dp
1929
1930 A_Row = 0.0_dp
1931 A_Col = 0.0_dp
1932
1933 f_Row = 0.0_dp
1934 f_Col = 0.0_dp
1935 f_Temp = 0.0_dp
1936
1937 t_Row = 0.0_dp
1938 t_Col = 0.0_dp
1939 t_Temp = 0.0_dp
1940
1941 pot_Row = 0.0_dp
1942 pot_Col = 0.0_dp
1943 pot_Temp = 0.0_dp
1944 ppot_Temp = 0.0_dp
1945
1946 frho_row = 0.0_dp
1947 frho_col = 0.0_dp
1948 rho_row = 0.0_dp
1949 rho_col = 0.0_dp
1950 rho_tmp = 0.0_dp
1951 dfrhodrho_row = 0.0_dp
1952 dfrhodrho_col = 0.0_dp
1953
1954 #endif
1955 rho = 0.0_dp
1956 frho = 0.0_dp
1957 dfrhodrho = 0.0_dp
1958
1959 end subroutine zero_work_arrays
1960
1961 function skipThisPair(atom1, atom2) result(skip_it)
1962 integer, intent(in) :: atom1
1963 integer, intent(in), optional :: atom2
1964 logical :: skip_it
1965 integer :: unique_id_1, unique_id_2
1966 integer :: me_i,me_j
1967 integer :: i
1968
1969 skip_it = .false.
1970
1971 !! there are a number of reasons to skip a pair or a particle
1972 !! mostly we do this to exclude atoms who are involved in short
1973 !! range interactions (bonds, bends, torsions), but we also need
1974 !! to exclude some overcounted interactions that result from
1975 !! the parallel decomposition
1976
1977 #ifdef IS_MPI
1978 !! in MPI, we have to look up the unique IDs for each atom
1979 unique_id_1 = AtomRowToGlobal(atom1)
1980 unique_id_2 = AtomColToGlobal(atom2)
1981 !! this situation should only arise in MPI simulations
1982 if (unique_id_1 == unique_id_2) then
1983 skip_it = .true.
1984 return
1985 end if
1986
1987 !! this prevents us from doing the pair on multiple processors
1988 if (unique_id_1 < unique_id_2) then
1989 if (mod(unique_id_1 + unique_id_2,2) == 0) then
1990 skip_it = .true.
1991 return
1992 endif
1993 else
1994 if (mod(unique_id_1 + unique_id_2,2) == 1) then
1995 skip_it = .true.
1996 return
1997 endif
1998 endif
1999 #else
2000 !! in the normal loop, the atom numbers are unique
2001 unique_id_1 = atom1
2002 unique_id_2 = atom2
2003 #endif
2004
2005 #ifdef IS_MPI
2006 do i = 1, nSkipsForRowAtom(atom1)
2007 if (skipsForRowAtom(atom1, i) .eq. unique_id_2) then
2008 skip_it = .true.
2009 return
2010 endif
2011 end do
2012 #else
2013 do i = 1, nSkipsForLocalAtom(atom1)
2014 if (skipsForLocalAtom(atom1, i) .eq. unique_id_2) then
2015 skip_it = .true.
2016 return
2017 endif
2018 end do
2019 #endif
2020
2021 return
2022 end function skipThisPair
2023
2024 function getTopoDistance(atom1, atom2) result(topoDist)
2025 integer, intent(in) :: atom1
2026 integer, intent(in) :: atom2
2027 integer :: topoDist
2028 integer :: unique_id_2
2029 integer :: i
2030
2031 #ifdef IS_MPI
2032 unique_id_2 = AtomColToGlobal(atom2)
2033 #else
2034 unique_id_2 = atom2
2035 #endif
2036
2037 ! zero is default for unconnected (i.e. normal) pair interactions
2038
2039 topoDist = 0
2040
2041 do i = 1, nTopoPairsForAtom(atom1)
2042 if (toposForAtom(atom1, i) .eq. unique_id_2) then
2043 topoDist = topoDistance(atom1, i)
2044 return
2045 endif
2046 end do
2047
2048 return
2049 end function getTopoDistance
2050
2051 function FF_UsesDirectionalAtoms() result(doesit)
2052 logical :: doesit
2053 doesit = FF_uses_DirectionalAtoms
2054 end function FF_UsesDirectionalAtoms
2055
2056 function FF_RequiresPrepairCalc() result(doesit)
2057 logical :: doesit
2058 doesit = FF_uses_EAM .or. FF_uses_SC
2059 end function FF_RequiresPrepairCalc
2060
2061 #ifdef PROFILE
2062 function getforcetime() result(totalforcetime)
2063 real(kind=dp) :: totalforcetime
2064 totalforcetime = forcetime
2065 end function getforcetime
2066 #endif
2067
2068 !! This cleans componets of force arrays belonging only to fortran
2069
2070 subroutine add_stress_tensor(dpair, fpair, tau)
2071
2072 real( kind = dp ), dimension(3), intent(in) :: dpair, fpair
2073 real( kind = dp ), dimension(9), intent(inout) :: tau
2074
2075 ! because the d vector is the rj - ri vector, and
2076 ! because fx, fy, fz are the force on atom i, we need a
2077 ! negative sign here:
2078
2079 tau(1) = tau(1) - dpair(1) * fpair(1)
2080 tau(2) = tau(2) - dpair(1) * fpair(2)
2081 tau(3) = tau(3) - dpair(1) * fpair(3)
2082 tau(4) = tau(4) - dpair(2) * fpair(1)
2083 tau(5) = tau(5) - dpair(2) * fpair(2)
2084 tau(6) = tau(6) - dpair(2) * fpair(3)
2085 tau(7) = tau(7) - dpair(3) * fpair(1)
2086 tau(8) = tau(8) - dpair(3) * fpair(2)
2087 tau(9) = tau(9) - dpair(3) * fpair(3)
2088
2089 end subroutine add_stress_tensor
2090
2091 end module doForces

Properties

Name Value
svn:keywords Author Id Revision Date