ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/doForces.F90
Revision: 1469
Committed: Mon Jul 19 14:07:59 2010 UTC (14 years, 9 months ago) by gezelter
File size: 68141 byte(s)
Log Message:
attempts at c++->fortran linkage.  Still busticated.

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

Properties

Name Value
svn:keywords Author Id Revision Date