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

Properties

Name Value
svn:keywords Author Id Revision Date