ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/doForces.F90
Revision: 1309
Committed: Tue Oct 21 18:23:31 2008 UTC (16 years, 7 months ago) by gezelter
File size: 60070 byte(s)
Log Message:
many changes to keep track of particle potentials for both bonded and non-bonded interactions

File Contents

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