ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/doForces.F90
Revision: 1489
Committed: Tue Aug 10 18:34:59 2010 UTC (14 years, 8 months ago) by gezelter
File size: 67487 byte(s)
Log Message:
Migrating Sutton-Chen from Fortran over to C++

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

Properties

Name Value
svn:keywords Author Id Revision Date