ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/doForces.F90
Revision: 1507
Committed: Tue Oct 5 18:34:55 2010 UTC (14 years, 7 months ago) by chuckv
File size: 67708 byte(s)
Log Message:
Replaced anint function in get_interatomic_vector with floor and ceiling functions.

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date