ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/doForces.F90
Revision: 1313
Committed: Wed Oct 22 20:01:49 2008 UTC (16 years, 6 months ago) by gezelter
File size: 60266 byte(s)
Log Message:
General bug-fixes and other changes to make particle pots work with
the Helfand Energy correlation function

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 I, 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.98 2008-10-22 20:01:48 gezelter Exp $, $Date: 2008-10-22 20:01:48 $, $Name: not supported by cvs2svn $, $Revision: 1.98 $
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 biggestAtypeCutoff = 0.0_dp
312 do i = 1, nAtypes
313 if (SimHasAtype(i)) then
314 call getElementProperty(atypes, i, "is_LennardJones", i_is_LJ)
315 call getElementProperty(atypes, i, "is_Electrostatic", i_is_Elect)
316 call getElementProperty(atypes, i, "is_Sticky", i_is_Sticky)
317 call getElementProperty(atypes, i, "is_StickyPower", i_is_StickyP)
318 call getElementProperty(atypes, i, "is_GayBerne", i_is_GB)
319 call getElementProperty(atypes, i, "is_EAM", i_is_EAM)
320 call getElementProperty(atypes, i, "is_Shape", i_is_Shape)
321 call getElementProperty(atypes, i, "is_SC", i_is_SC)
322
323 if (haveDefaultCutoffs) then
324 atypeMaxCutoff(i) = defaultRcut
325 else
326 if (i_is_LJ) then
327 thisRcut = getSigma(i) * 2.5_dp
328 if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
329 endif
330 if (i_is_Elect) then
331 thisRcut = defaultRcut
332 if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
333 endif
334 if (i_is_Sticky) then
335 thisRcut = getStickyCut(i)
336 if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
337 endif
338 if (i_is_StickyP) then
339 thisRcut = getStickyPowerCut(i)
340 if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
341 endif
342 if (i_is_GB) then
343 thisRcut = getGayBerneCut(i)
344 if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
345 endif
346 if (i_is_EAM) then
347 thisRcut = getEAMCut(i)
348 if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
349 endif
350 if (i_is_Shape) then
351 thisRcut = getShapeCut(i)
352 if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
353 endif
354 if (i_is_SC) then
355 thisRcut = getSCCut(i)
356 if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
357 endif
358 endif
359
360 if (atypeMaxCutoff(i).gt.biggestAtypeCutoff) then
361 biggestAtypeCutoff = atypeMaxCutoff(i)
362 endif
363
364 endif
365 enddo
366
367 istart = 1
368 jstart = 1
369 #ifdef IS_MPI
370 iend = nGroupsInRow
371 jend = nGroupsInCol
372 #else
373 iend = nGroups
374 jend = nGroups
375 #endif
376
377 !! allocate the groupToGtype and gtypeMaxCutoff here.
378 if(.not.allocated(groupToGtypeRow)) then
379 ! allocate(groupToGtype(iend))
380 allocate(groupToGtypeRow(iend))
381 else
382 deallocate(groupToGtypeRow)
383 allocate(groupToGtypeRow(iend))
384 endif
385 if(.not.allocated(groupMaxCutoffRow)) then
386 allocate(groupMaxCutoffRow(iend))
387 else
388 deallocate(groupMaxCutoffRow)
389 allocate(groupMaxCutoffRow(iend))
390 end if
391
392 if(.not.allocated(gtypeMaxCutoffRow)) then
393 allocate(gtypeMaxCutoffRow(iend))
394 else
395 deallocate(gtypeMaxCutoffRow)
396 allocate(gtypeMaxCutoffRow(iend))
397 endif
398
399
400 #ifdef IS_MPI
401 ! We only allocate new storage if we are in MPI because Ncol /= Nrow
402 if(.not.associated(groupToGtypeCol)) then
403 allocate(groupToGtypeCol(jend))
404 else
405 deallocate(groupToGtypeCol)
406 allocate(groupToGtypeCol(jend))
407 end if
408
409 if(.not.associated(groupMaxCutoffCol)) then
410 allocate(groupMaxCutoffCol(jend))
411 else
412 deallocate(groupMaxCutoffCol)
413 allocate(groupMaxCutoffCol(jend))
414 end if
415 if(.not.associated(gtypeMaxCutoffCol)) then
416 allocate(gtypeMaxCutoffCol(jend))
417 else
418 deallocate(gtypeMaxCutoffCol)
419 allocate(gtypeMaxCutoffCol(jend))
420 end if
421
422 groupMaxCutoffCol = 0.0_dp
423 gtypeMaxCutoffCol = 0.0_dp
424
425 #endif
426 groupMaxCutoffRow = 0.0_dp
427 gtypeMaxCutoffRow = 0.0_dp
428
429
430 !! first we do a single loop over the cutoff groups to find the
431 !! largest cutoff for any atypes present in this group. We also
432 !! create gtypes at this point.
433
434 tol = 1.0e-6_dp
435 nGroupTypesRow = 0
436 nGroupTypesCol = 0
437 do i = istart, iend
438 n_in_i = groupStartRow(i+1) - groupStartRow(i)
439 groupMaxCutoffRow(i) = 0.0_dp
440 do ia = groupStartRow(i), groupStartRow(i+1)-1
441 atom1 = groupListRow(ia)
442 #ifdef IS_MPI
443 me_i = atid_row(atom1)
444 #else
445 me_i = atid(atom1)
446 #endif
447 if (atypeMaxCutoff(me_i).gt.groupMaxCutoffRow(i)) then
448 groupMaxCutoffRow(i)=atypeMaxCutoff(me_i)
449 endif
450 enddo
451 if (nGroupTypesRow.eq.0) then
452 nGroupTypesRow = nGroupTypesRow + 1
453 gtypeMaxCutoffRow(nGroupTypesRow) = groupMaxCutoffRow(i)
454 groupToGtypeRow(i) = nGroupTypesRow
455 else
456 GtypeFound = .false.
457 do g = 1, nGroupTypesRow
458 if ( abs(groupMaxCutoffRow(i) - gtypeMaxCutoffRow(g)).lt.tol) then
459 groupToGtypeRow(i) = g
460 GtypeFound = .true.
461 endif
462 enddo
463 if (.not.GtypeFound) then
464 nGroupTypesRow = nGroupTypesRow + 1
465 gtypeMaxCutoffRow(nGroupTypesRow) = groupMaxCutoffRow(i)
466 groupToGtypeRow(i) = nGroupTypesRow
467 endif
468 endif
469 enddo
470
471 #ifdef IS_MPI
472 do j = jstart, jend
473 n_in_j = groupStartCol(j+1) - groupStartCol(j)
474 groupMaxCutoffCol(j) = 0.0_dp
475 do ja = groupStartCol(j), groupStartCol(j+1)-1
476 atom1 = groupListCol(ja)
477
478 me_j = atid_col(atom1)
479
480 if (atypeMaxCutoff(me_j).gt.groupMaxCutoffCol(j)) then
481 groupMaxCutoffCol(j)=atypeMaxCutoff(me_j)
482 endif
483 enddo
484
485 if (nGroupTypesCol.eq.0) then
486 nGroupTypesCol = nGroupTypesCol + 1
487 gtypeMaxCutoffCol(nGroupTypesCol) = groupMaxCutoffCol(j)
488 groupToGtypeCol(j) = nGroupTypesCol
489 else
490 GtypeFound = .false.
491 do g = 1, nGroupTypesCol
492 if ( abs(groupMaxCutoffCol(j) - gtypeMaxCutoffCol(g)).lt.tol) then
493 groupToGtypeCol(j) = g
494 GtypeFound = .true.
495 endif
496 enddo
497 if (.not.GtypeFound) then
498 nGroupTypesCol = nGroupTypesCol + 1
499 gtypeMaxCutoffCol(nGroupTypesCol) = groupMaxCutoffCol(j)
500 groupToGtypeCol(j) = nGroupTypesCol
501 endif
502 endif
503 enddo
504
505 #else
506 ! Set pointers to information we just found
507 nGroupTypesCol = nGroupTypesRow
508 groupToGtypeCol => groupToGtypeRow
509 gtypeMaxCutoffCol => gtypeMaxCutoffRow
510 groupMaxCutoffCol => groupMaxCutoffRow
511 #endif
512
513 !! allocate the gtypeCutoffMap here.
514 allocate(gtypeCutoffMap(nGroupTypesRow,nGroupTypesCol))
515 !! then we do a double loop over all the group TYPES to find the cutoff
516 !! map between groups of two types
517 tradRcut = max(maxval(gtypeMaxCutoffRow),maxval(gtypeMaxCutoffCol))
518
519 do i = 1, nGroupTypesRow
520 do j = 1, nGroupTypesCol
521
522 select case(cutoffPolicy)
523 case(TRADITIONAL_CUTOFF_POLICY)
524 thisRcut = tradRcut
525 case(MIX_CUTOFF_POLICY)
526 thisRcut = 0.5_dp * (gtypeMaxCutoffRow(i) + gtypeMaxCutoffCol(j))
527 case(MAX_CUTOFF_POLICY)
528 thisRcut = max(gtypeMaxCutoffRow(i), gtypeMaxCutoffCol(j))
529 case default
530 call handleError("createGtypeCutoffMap", "Unknown Cutoff Policy")
531 return
532 end select
533 gtypeCutoffMap(i,j)%rcut = thisRcut
534
535 if (thisRcut.gt.largestRcut) largestRcut = thisRcut
536
537 gtypeCutoffMap(i,j)%rcutsq = thisRcut*thisRcut
538
539 if (.not.haveSkinThickness) then
540 skinThickness = 1.0_dp
541 endif
542
543 gtypeCutoffMap(i,j)%rlistsq = (thisRcut + skinThickness)**2
544
545 ! sanity check
546
547 if (haveDefaultCutoffs) then
548 if (abs(gtypeCutoffMap(i,j)%rcut - defaultRcut).gt.0.0001) then
549 call handleError("createGtypeCutoffMap", "user-specified rCut does not match computed group Cutoff")
550 endif
551 endif
552 enddo
553 enddo
554
555 if(allocated(gtypeMaxCutoffRow)) deallocate(gtypeMaxCutoffRow)
556 if(allocated(groupMaxCutoffRow)) deallocate(groupMaxCutoffRow)
557 if(allocated(atypeMaxCutoff)) deallocate(atypeMaxCutoff)
558 #ifdef IS_MPI
559 if(associated(groupMaxCutoffCol)) deallocate(groupMaxCutoffCol)
560 if(associated(gtypeMaxCutoffCol)) deallocate(gtypeMaxCutoffCol)
561 #endif
562 groupMaxCutoffCol => null()
563 gtypeMaxCutoffCol => null()
564
565 haveGtypeCutoffMap = .true.
566 end subroutine createGtypeCutoffMap
567
568 subroutine setCutoffs(defRcut, defRsw, defSP, defSF)
569
570 real(kind=dp),intent(in) :: defRcut, defRsw
571 logical, intent(in) :: defSP, defSF
572 character(len = statusMsgSize) :: errMsg
573 integer :: localError
574
575 defaultRcut = defRcut
576 defaultRsw = defRsw
577
578 defaultDoShiftPot = defSP
579 defaultDoShiftFrc = defSF
580
581 if (abs(defaultRcut-defaultRsw) .lt. 0.0001) then
582 if (defaultDoShiftFrc) then
583 write(errMsg, *) &
584 'cutoffRadius and switchingRadius are set to the', newline &
585 // tab, 'same value. OOPSE will use shifted force', newline &
586 // tab, 'potentials instead of switching functions.'
587
588 call handleInfo("setCutoffs", errMsg)
589 else
590 write(errMsg, *) &
591 'cutoffRadius and switchingRadius are set to the', newline &
592 // tab, 'same value. OOPSE will use shifted', newline &
593 // tab, 'potentials instead of switching functions.'
594
595 call handleInfo("setCutoffs", errMsg)
596
597 defaultDoShiftPot = .true.
598 endif
599
600 endif
601
602 localError = 0
603 call setLJDefaultCutoff( defaultRcut, defaultDoShiftPot, &
604 defaultDoShiftFrc )
605 call setElectrostaticCutoffRadius( defaultRcut, defaultRsw )
606 call setCutoffEAM( defaultRcut )
607 call setCutoffSC( defaultRcut )
608 call setMnMDefaultCutoff( defaultRcut, defaultDoShiftPot, &
609 defaultDoShiftFrc )
610 call set_switch(defaultRsw, defaultRcut)
611 call setHmatDangerousRcutValue(defaultRcut)
612
613 haveDefaultCutoffs = .true.
614 haveGtypeCutoffMap = .false.
615
616 end subroutine setCutoffs
617
618 subroutine cWasLame()
619
620 VisitCutoffsAfterComputing = .true.
621 return
622
623 end subroutine cWasLame
624
625 subroutine setCutoffPolicy(cutPolicy)
626
627 integer, intent(in) :: cutPolicy
628
629 cutoffPolicy = cutPolicy
630 haveCutoffPolicy = .true.
631 haveGtypeCutoffMap = .false.
632
633 end subroutine setCutoffPolicy
634
635 subroutine setBoxDipole()
636
637 do_box_dipole = .true.
638
639 end subroutine setBoxDipole
640
641 subroutine getBoxDipole( box_dipole )
642
643 real(kind=dp), intent(inout), dimension(3) :: box_dipole
644
645 box_dipole = boxDipole
646
647 end subroutine getBoxDipole
648
649 subroutine setElectrostaticMethod( thisESM )
650
651 integer, intent(in) :: thisESM
652
653 electrostaticSummationMethod = thisESM
654 haveElectrostaticSummationMethod = .true.
655
656 end subroutine setElectrostaticMethod
657
658 subroutine setSkinThickness( thisSkin )
659
660 real(kind=dp), intent(in) :: thisSkin
661
662 skinThickness = thisSkin
663 haveSkinThickness = .true.
664 haveGtypeCutoffMap = .false.
665
666 end subroutine setSkinThickness
667
668 subroutine setSimVariables()
669 SIM_uses_DirectionalAtoms = SimUsesDirectionalAtoms()
670 SIM_uses_EAM = SimUsesEAM()
671 SIM_requires_postpair_calc = SimRequiresPostpairCalc()
672 SIM_requires_prepair_calc = SimRequiresPrepairCalc()
673 SIM_uses_PBC = SimUsesPBC()
674 SIM_uses_SC = SimUsesSC()
675 SIM_uses_AtomicVirial = SimUsesAtomicVirial()
676
677 haveSIMvariables = .true.
678
679 return
680 end subroutine setSimVariables
681
682 subroutine doReadyCheck(error)
683 integer, intent(out) :: error
684 integer :: myStatus
685
686 error = 0
687
688 if (.not. haveInteractionHash) then
689 call createInteractionHash()
690 endif
691
692 if (.not. haveGtypeCutoffMap) then
693 call createGtypeCutoffMap()
694 endif
695
696 if (VisitCutoffsAfterComputing) then
697 call set_switch(largestRcut, largestRcut)
698 call setHmatDangerousRcutValue(largestRcut)
699 call setCutoffEAM(largestRcut)
700 call setCutoffSC(largestRcut)
701 VisitCutoffsAfterComputing = .false.
702 endif
703
704 if (.not. haveSIMvariables) then
705 call setSimVariables()
706 endif
707
708 if (.not. haveNeighborList) then
709 write(default_error, *) 'neighbor list has not been initialized in doForces!'
710 error = -1
711 return
712 end if
713
714 if (.not. haveSaneForceField) then
715 write(default_error, *) 'Force Field is not sane in doForces!'
716 error = -1
717 return
718 end if
719
720 #ifdef IS_MPI
721 if (.not. isMPISimSet()) then
722 write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
723 error = -1
724 return
725 endif
726 #endif
727 return
728 end subroutine doReadyCheck
729
730
731 subroutine init_FF(thisStat)
732
733 integer, intent(out) :: thisStat
734 integer :: my_status, nMatches
735 integer, pointer :: MatchList(:) => null()
736
737 !! assume things are copacetic, unless they aren't
738 thisStat = 0
739
740 !! init_FF is called *after* all of the atom types have been
741 !! defined in atype_module using the new_atype subroutine.
742 !!
743 !! this will scan through the known atypes and figure out what
744 !! interactions are used by the force field.
745
746 FF_uses_DirectionalAtoms = .false.
747 FF_uses_Dipoles = .false.
748 FF_uses_GayBerne = .false.
749 FF_uses_EAM = .false.
750 FF_uses_SC = .false.
751
752 call getMatchingElementList(atypes, "is_Directional", .true., &
753 nMatches, MatchList)
754 if (nMatches .gt. 0) FF_uses_DirectionalAtoms = .true.
755
756 call getMatchingElementList(atypes, "is_Dipole", .true., &
757 nMatches, MatchList)
758 if (nMatches .gt. 0) FF_uses_Dipoles = .true.
759
760 call getMatchingElementList(atypes, "is_GayBerne", .true., &
761 nMatches, MatchList)
762 if (nMatches .gt. 0) FF_uses_GayBerne = .true.
763
764 call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList)
765 if (nMatches .gt. 0) FF_uses_EAM = .true.
766
767 call getMatchingElementList(atypes, "is_SC", .true., nMatches, MatchList)
768 if (nMatches .gt. 0) FF_uses_SC = .true.
769
770
771 haveSaneForceField = .true.
772
773 if (FF_uses_EAM) then
774 call init_EAM_FF(my_status)
775 if (my_status /= 0) then
776 write(default_error, *) "init_EAM_FF returned a bad status"
777 thisStat = -1
778 haveSaneForceField = .false.
779 return
780 end if
781 endif
782
783 if (.not. haveNeighborList) then
784 !! Create neighbor lists
785 call expandNeighborList(nLocal, my_status)
786 if (my_Status /= 0) then
787 write(default_error,*) "SimSetup: ExpandNeighborList returned error."
788 thisStat = -1
789 return
790 endif
791 haveNeighborList = .true.
792 endif
793
794 end subroutine init_FF
795
796
797 !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
798 !------------------------------------------------------------->
799 subroutine do_force_loop(q, q_group, A, eFrame, f, t, tau, pot, particle_pot, &
800 do_pot_c, do_stress_c, error)
801 !! Position array provided by C, dimensioned by getNlocal
802 real ( kind = dp ), dimension(3, nLocal) :: q
803 !! molecular center-of-mass position array
804 real ( kind = dp ), dimension(3, nGroups) :: q_group
805 !! Rotation Matrix for each long range particle in simulation.
806 real( kind = dp), dimension(9, nLocal) :: A
807 !! Unit vectors for dipoles (lab frame)
808 real( kind = dp ), dimension(9,nLocal) :: eFrame
809 !! Force array provided by C, dimensioned by getNlocal
810 real ( kind = dp ), dimension(3,nLocal) :: f
811 !! Torsion array provided by C, dimensioned by getNlocal
812 real( kind = dp ), dimension(3,nLocal) :: t
813
814 !! Stress Tensor
815 real( kind = dp), dimension(9) :: tau
816 real ( kind = dp ),dimension(LR_POT_TYPES) :: pot
817 real( kind = dp ), dimension(nLocal) :: particle_pot
818 logical ( kind = 2) :: do_pot_c, do_stress_c
819 logical :: do_pot
820 logical :: do_stress
821 logical :: in_switching_region
822 #ifdef IS_MPI
823 real( kind = DP ), dimension(LR_POT_TYPES) :: pot_local
824 integer :: nAtomsInRow
825 integer :: nAtomsInCol
826 integer :: nprocs
827 integer :: nGroupsInRow
828 integer :: nGroupsInCol
829 #endif
830 integer :: natoms
831 logical :: update_nlist
832 integer :: i, j, jstart, jend, jnab
833 integer :: istart, iend
834 integer :: ia, jb, atom1, atom2
835 integer :: nlist
836 real( kind = DP ) :: ratmsq, rgrpsq, rgrp, rag, vpair, vij
837 real( kind = DP ) :: sw, dswdr, swderiv, mf
838 real( kind = DP ) :: rVal
839 real(kind=dp),dimension(3) :: d_atm, d_grp, fpair, fij, fg, dag
840 real(kind=dp) :: rfpot, mu_i
841 real(kind=dp):: rCut
842 integer :: me_i, me_j, n_in_i, n_in_j
843 logical :: is_dp_i
844 integer :: neighborListSize
845 integer :: listerror, error
846 integer :: localError
847 integer :: propPack_i, propPack_j
848 integer :: loopStart, loopEnd, loop
849 integer :: iHash
850 integer :: i1, topoDist
851
852 !! the variables for the box dipole moment
853 #ifdef IS_MPI
854 integer :: pChgCount_local
855 integer :: nChgCount_local
856 real(kind=dp) :: pChg_local
857 real(kind=dp) :: nChg_local
858 real(kind=dp), dimension(3) :: pChgPos_local
859 real(kind=dp), dimension(3) :: nChgPos_local
860 real(kind=dp), dimension(3) :: dipVec_local
861 #endif
862 integer :: pChgCount
863 integer :: nChgCount
864 real(kind=dp) :: pChg
865 real(kind=dp) :: nChg
866 real(kind=dp) :: chg_value
867 real(kind=dp), dimension(3) :: pChgPos
868 real(kind=dp), dimension(3) :: nChgPos
869 real(kind=dp), dimension(3) :: dipVec
870 real(kind=dp), dimension(3) :: chgVec
871
872 !! initialize box dipole variables
873 if (do_box_dipole) then
874 #ifdef IS_MPI
875 pChg_local = 0.0_dp
876 nChg_local = 0.0_dp
877 pChgCount_local = 0
878 nChgCount_local = 0
879 do i=1, 3
880 pChgPos_local = 0.0_dp
881 nChgPos_local = 0.0_dp
882 dipVec_local = 0.0_dp
883 enddo
884 #endif
885 pChg = 0.0_dp
886 nChg = 0.0_dp
887 pChgCount = 0
888 nChgCount = 0
889 chg_value = 0.0_dp
890
891 do i=1, 3
892 pChgPos(i) = 0.0_dp
893 nChgPos(i) = 0.0_dp
894 dipVec(i) = 0.0_dp
895 chgVec(i) = 0.0_dp
896 boxDipole(i) = 0.0_dp
897 enddo
898 endif
899
900 !! initialize local variables
901
902 #ifdef IS_MPI
903 pot_local = 0.0_dp
904 nAtomsInRow = getNatomsInRow(plan_atom_row)
905 nAtomsInCol = getNatomsInCol(plan_atom_col)
906 nGroupsInRow = getNgroupsInRow(plan_group_row)
907 nGroupsInCol = getNgroupsInCol(plan_group_col)
908 #else
909 natoms = nlocal
910 #endif
911
912 call doReadyCheck(localError)
913 if ( localError .ne. 0 ) then
914 call handleError("do_force_loop", "Not Initialized")
915 error = -1
916 return
917 end if
918 call zero_work_arrays()
919
920 do_pot = do_pot_c
921 do_stress = do_stress_c
922
923 ! Gather all information needed by all force loops:
924
925 #ifdef IS_MPI
926
927 call gather(q, q_Row, plan_atom_row_3d)
928 call gather(q, q_Col, plan_atom_col_3d)
929
930 call gather(q_group, q_group_Row, plan_group_row_3d)
931 call gather(q_group, q_group_Col, plan_group_col_3d)
932
933 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
934 call gather(eFrame, eFrame_Row, plan_atom_row_rotation)
935 call gather(eFrame, eFrame_Col, plan_atom_col_rotation)
936
937 call gather(A, A_Row, plan_atom_row_rotation)
938 call gather(A, A_Col, plan_atom_col_rotation)
939 endif
940
941 #endif
942
943 !! Begin force loop timing:
944 #ifdef PROFILE
945 call cpu_time(forceTimeInitial)
946 nloops = nloops + 1
947 #endif
948
949 loopEnd = PAIR_LOOP
950 if (FF_RequiresPrepairCalc() .and. SIM_requires_prepair_calc) then
951 loopStart = PREPAIR_LOOP
952 else
953 loopStart = PAIR_LOOP
954 endif
955
956 do loop = loopStart, loopEnd
957
958 ! See if we need to update neighbor lists
959 ! (but only on the first time through):
960 if (loop .eq. loopStart) then
961 #ifdef IS_MPI
962 call checkNeighborList(nGroupsInRow, q_group_row, skinThickness, &
963 update_nlist)
964 #else
965 call checkNeighborList(nGroups, q_group, skinThickness, &
966 update_nlist)
967 #endif
968 endif
969
970 if (update_nlist) then
971 !! save current configuration and construct neighbor list
972 #ifdef IS_MPI
973 call saveNeighborList(nGroupsInRow, q_group_row)
974 #else
975 call saveNeighborList(nGroups, q_group)
976 #endif
977 neighborListSize = size(list)
978 nlist = 0
979 endif
980
981 istart = 1
982 #ifdef IS_MPI
983 iend = nGroupsInRow
984 #else
985 iend = nGroups - 1
986 #endif
987 outer: do i = istart, iend
988
989 if (update_nlist) point(i) = nlist + 1
990
991 n_in_i = groupStartRow(i+1) - groupStartRow(i)
992
993 if (update_nlist) then
994 #ifdef IS_MPI
995 jstart = 1
996 jend = nGroupsInCol
997 #else
998 jstart = i+1
999 jend = nGroups
1000 #endif
1001 else
1002 jstart = point(i)
1003 jend = point(i+1) - 1
1004 ! make sure group i has neighbors
1005 if (jstart .gt. jend) cycle outer
1006 endif
1007
1008 do jnab = jstart, jend
1009 if (update_nlist) then
1010 j = jnab
1011 else
1012 j = list(jnab)
1013 endif
1014
1015 #ifdef IS_MPI
1016 me_j = atid_col(j)
1017 call get_interatomic_vector(q_group_Row(:,i), &
1018 q_group_Col(:,j), d_grp, rgrpsq)
1019 #else
1020 me_j = atid(j)
1021 call get_interatomic_vector(q_group(:,i), &
1022 q_group(:,j), d_grp, rgrpsq)
1023 #endif
1024
1025 if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rListsq) then
1026 if (update_nlist) then
1027 nlist = nlist + 1
1028
1029 if (nlist > neighborListSize) then
1030 #ifdef IS_MPI
1031 call expandNeighborList(nGroupsInRow, listerror)
1032 #else
1033 call expandNeighborList(nGroups, listerror)
1034 #endif
1035 if (listerror /= 0) then
1036 error = -1
1037 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
1038 return
1039 end if
1040 neighborListSize = size(list)
1041 endif
1042
1043 list(nlist) = j
1044 endif
1045
1046 if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCutsq) then
1047
1048 rCut = gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCut
1049 if (loop .eq. PAIR_LOOP) then
1050 vij = 0.0_dp
1051 fij(1) = 0.0_dp
1052 fij(2) = 0.0_dp
1053 fij(3) = 0.0_dp
1054 endif
1055
1056 call get_switch(rgrpsq, sw, dswdr,rgrp, in_switching_region)
1057
1058 n_in_j = groupStartCol(j+1) - groupStartCol(j)
1059
1060 do ia = groupStartRow(i), groupStartRow(i+1)-1
1061
1062 atom1 = groupListRow(ia)
1063
1064 inner: do jb = groupStartCol(j), groupStartCol(j+1)-1
1065
1066 atom2 = groupListCol(jb)
1067
1068 if (skipThisPair(atom1, atom2)) cycle inner
1069
1070 if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
1071 d_atm(1) = d_grp(1)
1072 d_atm(2) = d_grp(2)
1073 d_atm(3) = d_grp(3)
1074 ratmsq = rgrpsq
1075 else
1076 #ifdef IS_MPI
1077 call get_interatomic_vector(q_Row(:,atom1), &
1078 q_Col(:,atom2), d_atm, ratmsq)
1079 #else
1080 call get_interatomic_vector(q(:,atom1), &
1081 q(:,atom2), d_atm, ratmsq)
1082 #endif
1083 endif
1084
1085 topoDist = getTopoDistance(atom1, atom2)
1086
1087 if (loop .eq. PREPAIR_LOOP) then
1088 #ifdef IS_MPI
1089 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
1090 rgrpsq, d_grp, rCut, do_pot, do_stress, &
1091 eFrame, A, f, t, pot_local)
1092 #else
1093 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
1094 rgrpsq, d_grp, rCut, do_pot, do_stress, &
1095 eFrame, A, f, t, pot)
1096 #endif
1097 else
1098 #ifdef IS_MPI
1099 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
1100 do_pot, eFrame, A, f, t, pot_local, particle_pot, vpair, &
1101 fpair, d_grp, rgrp, rCut, topoDist)
1102 ! particle_pot will be accumulated from row & column
1103 ! arrays later
1104 #else
1105 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
1106 do_pot, eFrame, A, f, t, pot, particle_pot, vpair, &
1107 fpair, d_grp, rgrp, rCut, topoDist)
1108 #endif
1109 vij = vij + vpair
1110 fij(1) = fij(1) + fpair(1)
1111 fij(2) = fij(2) + fpair(2)
1112 fij(3) = fij(3) + fpair(3)
1113 if (do_stress) then
1114 call add_stress_tensor(d_atm, fpair, tau)
1115 endif
1116 endif
1117 enddo inner
1118 enddo
1119
1120 if (loop .eq. PAIR_LOOP) then
1121 if (in_switching_region) then
1122 swderiv = vij*dswdr/rgrp
1123 fg = swderiv*d_grp
1124
1125 fij(1) = fij(1) + fg(1)
1126 fij(2) = fij(2) + fg(2)
1127 fij(3) = fij(3) + fg(3)
1128
1129 if (do_stress .and. (n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
1130 call add_stress_tensor(d_atm, fg, tau)
1131 endif
1132
1133 do ia=groupStartRow(i), groupStartRow(i+1)-1
1134 atom1=groupListRow(ia)
1135 mf = mfactRow(atom1)
1136 ! fg is the force on atom ia due to cutoff group's
1137 ! presence in switching region
1138 fg = swderiv*d_grp*mf
1139 #ifdef IS_MPI
1140 f_Row(1,atom1) = f_Row(1,atom1) + fg(1)
1141 f_Row(2,atom1) = f_Row(2,atom1) + fg(2)
1142 f_Row(3,atom1) = f_Row(3,atom1) + fg(3)
1143 #else
1144 f(1,atom1) = f(1,atom1) + fg(1)
1145 f(2,atom1) = f(2,atom1) + fg(2)
1146 f(3,atom1) = f(3,atom1) + fg(3)
1147 #endif
1148 if (n_in_i .gt. 1) then
1149 if (do_stress.and.SIM_uses_AtomicVirial) then
1150 ! find the distance between the atom and the center of
1151 ! the cutoff group:
1152 #ifdef IS_MPI
1153 call get_interatomic_vector(q_Row(:,atom1), &
1154 q_group_Row(:,i), dag, rag)
1155 #else
1156 call get_interatomic_vector(q(:,atom1), &
1157 q_group(:,i), dag, rag)
1158 #endif
1159 call add_stress_tensor(dag,fg,tau)
1160 endif
1161 endif
1162 enddo
1163
1164 do jb=groupStartCol(j), groupStartCol(j+1)-1
1165 atom2=groupListCol(jb)
1166 mf = mfactCol(atom2)
1167 ! fg is the force on atom jb due to cutoff group's
1168 ! presence in switching region
1169 fg = -swderiv*d_grp*mf
1170 #ifdef IS_MPI
1171 f_Col(1,atom2) = f_Col(1,atom2) + fg(1)
1172 f_Col(2,atom2) = f_Col(2,atom2) + fg(2)
1173 f_Col(3,atom2) = f_Col(3,atom2) + fg(3)
1174 #else
1175 f(1,atom2) = f(1,atom2) + fg(1)
1176 f(2,atom2) = f(2,atom2) + fg(2)
1177 f(3,atom2) = f(3,atom2) + fg(3)
1178 #endif
1179 if (n_in_j .gt. 1) then
1180 if (do_stress.and.SIM_uses_AtomicVirial) then
1181 ! find the distance between the atom and the center of
1182 ! the cutoff group:
1183 #ifdef IS_MPI
1184 call get_interatomic_vector(q_Col(:,atom2), &
1185 q_group_Col(:,j), dag, rag)
1186 #else
1187 call get_interatomic_vector(q(:,atom2), &
1188 q_group(:,j), dag, rag)
1189 #endif
1190 call add_stress_tensor(dag,fg,tau)
1191 endif
1192 endif
1193 enddo
1194 endif
1195 !if (do_stress.and.(.not.SIM_uses_AtomicVirial)) then
1196 ! call add_stress_tensor(d_grp, fij, tau)
1197 !endif
1198 endif
1199 endif
1200 endif
1201 enddo
1202
1203 enddo outer
1204
1205 if (update_nlist) then
1206 #ifdef IS_MPI
1207 point(nGroupsInRow + 1) = nlist + 1
1208 #else
1209 point(nGroups) = nlist + 1
1210 #endif
1211 if (loop .eq. PREPAIR_LOOP) then
1212 ! we just did the neighbor list update on the first
1213 ! pass, so we don't need to do it
1214 ! again on the second pass
1215 update_nlist = .false.
1216 endif
1217 endif
1218
1219 if (loop .eq. PREPAIR_LOOP) then
1220 #ifdef IS_MPI
1221 call do_preforce(nlocal, pot_local, particle_pot)
1222 #else
1223 call do_preforce(nlocal, pot, particle_pot)
1224 #endif
1225 endif
1226
1227 enddo
1228
1229 !! Do timing
1230 #ifdef PROFILE
1231 call cpu_time(forceTimeFinal)
1232 forceTime = forceTime + forceTimeFinal - forceTimeInitial
1233 #endif
1234
1235 #ifdef IS_MPI
1236 !!distribute forces
1237
1238 f_temp = 0.0_dp
1239 call scatter(f_Row,f_temp,plan_atom_row_3d)
1240 do i = 1,nlocal
1241 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1242 end do
1243
1244 f_temp = 0.0_dp
1245 call scatter(f_Col,f_temp,plan_atom_col_3d)
1246 do i = 1,nlocal
1247 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1248 end do
1249
1250 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
1251 t_temp = 0.0_dp
1252 call scatter(t_Row,t_temp,plan_atom_row_3d)
1253 do i = 1,nlocal
1254 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1255 end do
1256 t_temp = 0.0_dp
1257 call scatter(t_Col,t_temp,plan_atom_col_3d)
1258
1259 do i = 1,nlocal
1260 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1261 end do
1262 endif
1263
1264 if (do_pot) then
1265 ! scatter/gather pot_row into the members of my column
1266 do i = 1,LR_POT_TYPES
1267 call scatter(pot_Row(i,:), pot_Temp(i,:), plan_atom_row)
1268 end do
1269 ! scatter/gather pot_local into all other procs
1270 ! add resultant to get total pot
1271 do i = 1, nlocal
1272 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES) &
1273 + pot_Temp(1:LR_POT_TYPES,i)
1274 enddo
1275
1276 do i = 1,LR_POT_TYPES
1277 particle_pot(1:nlocal) = particle_pot(1:nlocal) + pot_Temp(i,1:nlocal)
1278 enddo
1279
1280 pot_Temp = 0.0_DP
1281
1282 do i = 1,LR_POT_TYPES
1283 call scatter(pot_Col(i,:), pot_Temp(i,:), plan_atom_col)
1284 end do
1285
1286 do i = 1, nlocal
1287 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES)&
1288 + pot_Temp(1:LR_POT_TYPES,i)
1289 enddo
1290
1291 do i = 1,LR_POT_TYPES
1292 particle_pot(1:nlocal) = particle_pot(1:nlocal) + pot_Temp(i,1:nlocal)
1293 enddo
1294
1295 ppot_Temp = 0.0_DP
1296
1297 call scatter(ppot_Row(:), ppot_Temp(:), plan_atom_row)
1298 do i = 1, nlocal
1299 particle_pot(i) = particle_pot(i) + ppot_Temp(i)
1300 enddo
1301
1302 ppot_Temp = 0.0_DP
1303
1304 call scatter(ppot_Col(:), ppot_Temp(:), plan_atom_col)
1305 do i = 1, nlocal
1306 particle_pot(i) = particle_pot(i) + ppot_Temp(i)
1307 enddo
1308
1309
1310 endif
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 #ifdef IS_MPI
1327 call self_self(i, eFrame, pot_local(ELECTROSTATIC_POT), &
1328 t, do_pot)
1329 #else
1330 call self_self(i, eFrame, pot(ELECTROSTATIC_POT), &
1331 t, do_pot)
1332 #endif
1333 endif
1334
1335
1336 if (electrostaticSummationMethod.eq.REACTION_FIELD) then
1337
1338 ! loop over the excludes to accumulate RF stuff we've
1339 ! left out of the normal pair loop
1340
1341 do i1 = 1, nSkipsForAtom(i)
1342 j = skipsForAtom(i, i1)
1343
1344 ! prevent overcounting of the skips
1345 if (i.lt.j) then
1346 call get_interatomic_vector(q(:,i), q(:,j), d_atm, ratmsq)
1347 rVal = sqrt(ratmsq)
1348 call get_switch(ratmsq, sw, dswdr, rVal,in_switching_region)
1349 #ifdef IS_MPI
1350 call rf_self_excludes(i, j, sw, 1.0_dp, eFrame, d_atm, rVal, &
1351 vpair, pot_local(ELECTROSTATIC_POT), f, t, do_pot)
1352 #else
1353 call rf_self_excludes(i, j, sw, 1.0_dp, eFrame, d_atm, rVal, &
1354 vpair, pot(ELECTROSTATIC_POT), f, t, do_pot)
1355 #endif
1356 endif
1357 enddo
1358 endif
1359
1360 if (do_box_dipole) then
1361 #ifdef IS_MPI
1362 call accumulate_box_dipole(i, eFrame, q(:,i), pChg_local, &
1363 nChg_local, pChgPos_local, nChgPos_local, dipVec_local, &
1364 pChgCount_local, nChgCount_local)
1365 #else
1366 call accumulate_box_dipole(i, eFrame, q(:,i), pChg, nChg, &
1367 pChgPos, nChgPos, dipVec, pChgCount, nChgCount)
1368 #endif
1369 endif
1370 enddo
1371 endif
1372
1373 #ifdef IS_MPI
1374 if (do_pot) then
1375 #ifdef SINGLE_PRECISION
1376 call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_real,mpi_sum, &
1377 mpi_comm_world,mpi_err)
1378 #else
1379 call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_double_precision, &
1380 mpi_sum, mpi_comm_world,mpi_err)
1381 #endif
1382 endif
1383
1384 if (do_box_dipole) then
1385
1386 #ifdef SINGLE_PRECISION
1387 call mpi_allreduce(pChg_local, pChg, 1, mpi_real, mpi_sum, &
1388 mpi_comm_world, mpi_err)
1389 call mpi_allreduce(nChg_local, nChg, 1, mpi_real, mpi_sum, &
1390 mpi_comm_world, mpi_err)
1391 call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer, mpi_sum,&
1392 mpi_comm_world, mpi_err)
1393 call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer, mpi_sum,&
1394 mpi_comm_world, mpi_err)
1395 call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_real, mpi_sum, &
1396 mpi_comm_world, mpi_err)
1397 call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_real, mpi_sum, &
1398 mpi_comm_world, mpi_err)
1399 call mpi_allreduce(dipVec_local, dipVec, 3, mpi_real, mpi_sum, &
1400 mpi_comm_world, mpi_err)
1401 #else
1402 call mpi_allreduce(pChg_local, pChg, 1, mpi_double_precision, mpi_sum, &
1403 mpi_comm_world, mpi_err)
1404 call mpi_allreduce(nChg_local, nChg, 1, mpi_double_precision, mpi_sum, &
1405 mpi_comm_world, mpi_err)
1406 call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer,&
1407 mpi_sum, mpi_comm_world, mpi_err)
1408 call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer,&
1409 mpi_sum, mpi_comm_world, mpi_err)
1410 call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_double_precision, &
1411 mpi_sum, mpi_comm_world, mpi_err)
1412 call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_double_precision, &
1413 mpi_sum, mpi_comm_world, mpi_err)
1414 call mpi_allreduce(dipVec_local, dipVec, 3, mpi_double_precision, &
1415 mpi_sum, mpi_comm_world, mpi_err)
1416 #endif
1417
1418 endif
1419
1420 #endif
1421
1422 if (do_box_dipole) then
1423 ! first load the accumulated dipole moment (if dipoles were present)
1424 boxDipole(1) = dipVec(1)
1425 boxDipole(2) = dipVec(2)
1426 boxDipole(3) = dipVec(3)
1427
1428 ! now include the dipole moment due to charges
1429 ! use the lesser of the positive and negative charge totals
1430 if (nChg .le. pChg) then
1431 chg_value = nChg
1432 else
1433 chg_value = pChg
1434 endif
1435
1436 ! find the average positions
1437 if (pChgCount .gt. 0 .and. nChgCount .gt. 0) then
1438 pChgPos = pChgPos / pChgCount
1439 nChgPos = nChgPos / nChgCount
1440 endif
1441
1442 ! dipole is from the negative to the positive (physics notation)
1443 chgVec(1) = pChgPos(1) - nChgPos(1)
1444 chgVec(2) = pChgPos(2) - nChgPos(2)
1445 chgVec(3) = pChgPos(3) - nChgPos(3)
1446
1447 boxDipole(1) = boxDipole(1) + chgVec(1) * chg_value
1448 boxDipole(2) = boxDipole(2) + chgVec(2) * chg_value
1449 boxDipole(3) = boxDipole(3) + chgVec(3) * chg_value
1450
1451 endif
1452
1453 end subroutine do_force_loop
1454
1455 subroutine do_pair(i, j, rijsq, d, sw, do_pot, &
1456 eFrame, A, f, t, pot, particle_pot, vpair, &
1457 fpair, d_grp, r_grp, rCut, topoDist)
1458
1459 real( kind = dp ) :: vpair, sw
1460 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1461 real( kind = dp ), dimension(nLocal) :: particle_pot
1462 real( kind = dp ), dimension(3) :: fpair
1463 real( kind = dp ), dimension(nLocal) :: mfact
1464 real( kind = dp ), dimension(9,nLocal) :: eFrame
1465 real( kind = dp ), dimension(9,nLocal) :: A
1466 real( kind = dp ), dimension(3,nLocal) :: f
1467 real( kind = dp ), dimension(3,nLocal) :: t
1468
1469 logical, intent(inout) :: do_pot
1470 integer, intent(in) :: i, j
1471 real ( kind = dp ), intent(inout) :: rijsq
1472 real ( kind = dp ), intent(inout) :: r_grp
1473 real ( kind = dp ), intent(inout) :: d(3)
1474 real ( kind = dp ), intent(inout) :: d_grp(3)
1475 real ( kind = dp ), intent(inout) :: rCut
1476 integer, intent(inout) :: topoDist
1477 real ( kind = dp ) :: r, pair_pot, vdwMult, electroMult
1478 real ( kind = dp ) :: a_k, b_k, c_k, d_k, dx
1479 integer :: me_i, me_j
1480 integer :: k
1481
1482 integer :: iHash
1483
1484 !!$ write(*,*) i, j, rijsq, d(1), d(2), d(3), sw, do_pot
1485 !!$ write(*,*) particle_pot(1), vpair, fpair(1), fpair(2), fpair(3)
1486 !!$ write(*,*) rCut
1487
1488
1489 r = sqrt(rijsq)
1490
1491 vpair = 0.0_dp
1492 fpair(1:3) = 0.0_dp
1493
1494 #ifdef IS_MPI
1495 me_i = atid_row(i)
1496 me_j = atid_col(j)
1497 #else
1498 me_i = atid(i)
1499 me_j = atid(j)
1500 #endif
1501
1502 iHash = InteractionHash(me_i, me_j)
1503
1504 vdwMult = vdwScale(topoDist)
1505 electroMult = electrostaticScale(topoDist)
1506
1507 if ( iand(iHash, LJ_PAIR).ne.0 ) then
1508 call do_lj_pair(i, j, d, r, rijsq, rcut, sw, vdwMult, vpair, fpair, &
1509 pot(VDW_POT), f, do_pot)
1510 endif
1511
1512 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1513 call doElectrostaticPair(i, j, d, r, rijsq, rcut, sw, electroMult, &
1514 vpair, fpair, pot(ELECTROSTATIC_POT), eFrame, f, t, do_pot)
1515 endif
1516
1517 if ( iand(iHash, STICKY_PAIR).ne.0 ) then
1518 call do_sticky_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1519 pot(HB_POT), A, f, t, do_pot)
1520 endif
1521
1522 if ( iand(iHash, STICKYPOWER_PAIR).ne.0 ) then
1523 call do_sticky_power_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1524 pot(HB_POT), A, f, t, do_pot)
1525 endif
1526
1527 if ( iand(iHash, GAYBERNE_PAIR).ne.0 ) then
1528 call do_gb_pair(i, j, d, r, rijsq, sw, vdwMult, vpair, fpair, &
1529 pot(VDW_POT), A, f, t, do_pot)
1530 endif
1531
1532 if ( iand(iHash, GAYBERNE_LJ).ne.0 ) then
1533 call do_gb_pair(i, j, d, r, rijsq, sw, vdwMult, vpair, fpair, &
1534 pot(VDW_POT), A, f, t, do_pot)
1535 endif
1536
1537 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1538 call do_eam_pair(i, j, d, r, rijsq, sw, vpair, particle_pot, &
1539 fpair, pot(METALLIC_POT), f, do_pot)
1540 endif
1541
1542 if ( iand(iHash, SHAPE_PAIR).ne.0 ) then
1543 call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1544 pot(VDW_POT), A, f, t, do_pot)
1545 endif
1546
1547 if ( iand(iHash, SHAPE_LJ).ne.0 ) then
1548 call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1549 pot(VDW_POT), A, f, t, do_pot)
1550 endif
1551
1552 if ( iand(iHash, SC_PAIR).ne.0 ) then
1553 call do_SC_pair(i, j, d, r, rijsq, rcut, sw, vpair, particle_pot, &
1554 fpair, pot(METALLIC_POT), f, do_pot)
1555 endif
1556
1557 if ( iand(iHash, MNM_PAIR).ne.0 ) then
1558 call do_mnm_pair(i, j, d, r, rijsq, rcut, sw, vdwMult, vpair, fpair, &
1559 pot(VDW_POT), A, f, t, do_pot)
1560 endif
1561 !!$
1562 !!$ particle_pot(i) = particle_pot(i) + vpair*sw
1563 !!$ particle_pot(j) = particle_pot(j) + vpair*sw
1564
1565 end subroutine do_pair
1566
1567 subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, rCut, &
1568 do_pot, do_stress, eFrame, A, f, t, pot)
1569
1570 real( kind = dp ) :: sw
1571 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1572 real( kind = dp ), dimension(9,nLocal) :: eFrame
1573 real (kind=dp), dimension(9,nLocal) :: A
1574 real (kind=dp), dimension(3,nLocal) :: f
1575 real (kind=dp), dimension(3,nLocal) :: t
1576
1577 logical, intent(inout) :: do_pot, do_stress
1578 integer, intent(in) :: i, j
1579 real ( kind = dp ), intent(inout) :: rijsq, rcijsq, rCut
1580 real ( kind = dp ) :: r, rc
1581 real ( kind = dp ), intent(inout) :: d(3), dc(3)
1582
1583 integer :: me_i, me_j, iHash
1584
1585 r = sqrt(rijsq)
1586
1587 #ifdef IS_MPI
1588 me_i = atid_row(i)
1589 me_j = atid_col(j)
1590 #else
1591 me_i = atid(i)
1592 me_j = atid(j)
1593 #endif
1594
1595 iHash = InteractionHash(me_i, me_j)
1596
1597 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1598 call calc_EAM_prepair_rho(i, j, d, r, rijsq)
1599 endif
1600
1601 if ( iand(iHash, SC_PAIR).ne.0 ) then
1602 call calc_SC_prepair_rho(i, j, d, r, rijsq, rcut )
1603 endif
1604
1605 end subroutine do_prepair
1606
1607
1608 subroutine do_preforce(nlocal, pot, particle_pot)
1609 integer :: nlocal
1610 real( kind = dp ),dimension(LR_POT_TYPES) :: pot
1611 real( kind = dp ),dimension(nlocal) :: particle_pot
1612
1613 if (FF_uses_EAM .and. SIM_uses_EAM) then
1614 call calc_EAM_preforce_Frho(nlocal, pot(METALLIC_POT), particle_pot)
1615 endif
1616 if (FF_uses_SC .and. SIM_uses_SC) then
1617 call calc_SC_preforce_Frho(nlocal, pot(METALLIC_POT), particle_pot)
1618 endif
1619 end subroutine do_preforce
1620
1621
1622 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
1623
1624 real (kind = dp), dimension(3) :: q_i
1625 real (kind = dp), dimension(3) :: q_j
1626 real ( kind = dp ), intent(out) :: r_sq
1627 real( kind = dp ) :: d(3), scaled(3)
1628 integer i
1629
1630 d(1) = q_j(1) - q_i(1)
1631 d(2) = q_j(2) - q_i(2)
1632 d(3) = q_j(3) - q_i(3)
1633
1634 ! Wrap back into periodic box if necessary
1635 if ( SIM_uses_PBC ) then
1636
1637 if( .not.boxIsOrthorhombic ) then
1638 ! calc the scaled coordinates.
1639 ! scaled = matmul(HmatInv, d)
1640
1641 scaled(1) = HmatInv(1,1)*d(1) + HmatInv(1,2)*d(2) + HmatInv(1,3)*d(3)
1642 scaled(2) = HmatInv(2,1)*d(1) + HmatInv(2,2)*d(2) + HmatInv(2,3)*d(3)
1643 scaled(3) = HmatInv(3,1)*d(1) + HmatInv(3,2)*d(2) + HmatInv(3,3)*d(3)
1644
1645 ! wrap the scaled coordinates
1646
1647 scaled(1) = scaled(1) - anint(scaled(1), kind=dp)
1648 scaled(2) = scaled(2) - anint(scaled(2), kind=dp)
1649 scaled(3) = scaled(3) - anint(scaled(3), kind=dp)
1650
1651 ! calc the wrapped real coordinates from the wrapped scaled
1652 ! coordinates
1653 ! d = matmul(Hmat,scaled)
1654 d(1)= Hmat(1,1)*scaled(1) + Hmat(1,2)*scaled(2) + Hmat(1,3)*scaled(3)
1655 d(2)= Hmat(2,1)*scaled(1) + Hmat(2,2)*scaled(2) + Hmat(2,3)*scaled(3)
1656 d(3)= Hmat(3,1)*scaled(1) + Hmat(3,2)*scaled(2) + Hmat(3,3)*scaled(3)
1657
1658 else
1659 ! calc the scaled coordinates.
1660
1661 scaled(1) = d(1) * HmatInv(1,1)
1662 scaled(2) = d(2) * HmatInv(2,2)
1663 scaled(3) = d(3) * HmatInv(3,3)
1664
1665 ! wrap the scaled coordinates
1666
1667 scaled(1) = scaled(1) - anint(scaled(1), kind=dp)
1668 scaled(2) = scaled(2) - anint(scaled(2), kind=dp)
1669 scaled(3) = scaled(3) - anint(scaled(3), kind=dp)
1670
1671 ! calc the wrapped real coordinates from the wrapped scaled
1672 ! coordinates
1673
1674 d(1) = scaled(1)*Hmat(1,1)
1675 d(2) = scaled(2)*Hmat(2,2)
1676 d(3) = scaled(3)*Hmat(3,3)
1677
1678 endif
1679
1680 endif
1681
1682 r_sq = d(1)*d(1) + d(2)*d(2) + d(3)*d(3)
1683
1684 end subroutine get_interatomic_vector
1685
1686 subroutine zero_work_arrays()
1687
1688 #ifdef IS_MPI
1689
1690 q_Row = 0.0_dp
1691 q_Col = 0.0_dp
1692
1693 q_group_Row = 0.0_dp
1694 q_group_Col = 0.0_dp
1695
1696 eFrame_Row = 0.0_dp
1697 eFrame_Col = 0.0_dp
1698
1699 A_Row = 0.0_dp
1700 A_Col = 0.0_dp
1701
1702 f_Row = 0.0_dp
1703 f_Col = 0.0_dp
1704 f_Temp = 0.0_dp
1705
1706 t_Row = 0.0_dp
1707 t_Col = 0.0_dp
1708 t_Temp = 0.0_dp
1709
1710 pot_Row = 0.0_dp
1711 pot_Col = 0.0_dp
1712 pot_Temp = 0.0_dp
1713 ppot_Temp = 0.0_dp
1714
1715 #endif
1716
1717 if (FF_uses_EAM .and. SIM_uses_EAM) then
1718 call clean_EAM()
1719 endif
1720
1721 end subroutine zero_work_arrays
1722
1723 function skipThisPair(atom1, atom2) result(skip_it)
1724 integer, intent(in) :: atom1
1725 integer, intent(in), optional :: atom2
1726 logical :: skip_it
1727 integer :: unique_id_1, unique_id_2
1728 integer :: me_i,me_j
1729 integer :: i
1730
1731 skip_it = .false.
1732
1733 !! there are a number of reasons to skip a pair or a particle
1734 !! mostly we do this to exclude atoms who are involved in short
1735 !! range interactions (bonds, bends, torsions), but we also need
1736 !! to exclude some overcounted interactions that result from
1737 !! the parallel decomposition
1738
1739 #ifdef IS_MPI
1740 !! in MPI, we have to look up the unique IDs for each atom
1741 unique_id_1 = AtomRowToGlobal(atom1)
1742 unique_id_2 = AtomColToGlobal(atom2)
1743 !! this situation should only arise in MPI simulations
1744 if (unique_id_1 == unique_id_2) then
1745 skip_it = .true.
1746 return
1747 end if
1748
1749 !! this prevents us from doing the pair on multiple processors
1750 if (unique_id_1 < unique_id_2) then
1751 if (mod(unique_id_1 + unique_id_2,2) == 0) then
1752 skip_it = .true.
1753 return
1754 endif
1755 else
1756 if (mod(unique_id_1 + unique_id_2,2) == 1) then
1757 skip_it = .true.
1758 return
1759 endif
1760 endif
1761 #else
1762 !! in the normal loop, the atom numbers are unique
1763 unique_id_1 = atom1
1764 unique_id_2 = atom2
1765 #endif
1766
1767 do i = 1, nSkipsForAtom(atom1)
1768 if (skipsForAtom(atom1, i) .eq. unique_id_2) then
1769 skip_it = .true.
1770 return
1771 endif
1772 end do
1773
1774 return
1775 end function skipThisPair
1776
1777 function getTopoDistance(atom1, atom2) result(topoDist)
1778 integer, intent(in) :: atom1
1779 integer, intent(in) :: atom2
1780 integer :: topoDist
1781 integer :: unique_id_2
1782 integer :: i
1783
1784 #ifdef IS_MPI
1785 unique_id_2 = AtomColToGlobal(atom2)
1786 #else
1787 unique_id_2 = atom2
1788 #endif
1789
1790 ! zero is default for unconnected (i.e. normal) pair interactions
1791
1792 topoDist = 0
1793
1794 do i = 1, nTopoPairsForAtom(atom1)
1795 if (toposForAtom(atom1, i) .eq. unique_id_2) then
1796 topoDist = topoDistance(atom1, i)
1797 return
1798 endif
1799 end do
1800
1801 return
1802 end function getTopoDistance
1803
1804 function FF_UsesDirectionalAtoms() result(doesit)
1805 logical :: doesit
1806 doesit = FF_uses_DirectionalAtoms
1807 end function FF_UsesDirectionalAtoms
1808
1809 function FF_RequiresPrepairCalc() result(doesit)
1810 logical :: doesit
1811 doesit = FF_uses_EAM .or. FF_uses_SC
1812 end function FF_RequiresPrepairCalc
1813
1814 #ifdef PROFILE
1815 function getforcetime() result(totalforcetime)
1816 real(kind=dp) :: totalforcetime
1817 totalforcetime = forcetime
1818 end function getforcetime
1819 #endif
1820
1821 !! This cleans componets of force arrays belonging only to fortran
1822
1823 subroutine add_stress_tensor(dpair, fpair, tau)
1824
1825 real( kind = dp ), dimension(3), intent(in) :: dpair, fpair
1826 real( kind = dp ), dimension(9), intent(inout) :: tau
1827
1828 ! because the d vector is the rj - ri vector, and
1829 ! because fx, fy, fz are the force on atom i, we need a
1830 ! negative sign here:
1831
1832 tau(1) = tau(1) - dpair(1) * fpair(1)
1833 tau(2) = tau(2) - dpair(1) * fpair(2)
1834 tau(3) = tau(3) - dpair(1) * fpair(3)
1835 tau(4) = tau(4) - dpair(2) * fpair(1)
1836 tau(5) = tau(5) - dpair(2) * fpair(2)
1837 tau(6) = tau(6) - dpair(2) * fpair(3)
1838 tau(7) = tau(7) - dpair(3) * fpair(1)
1839 tau(8) = tau(8) - dpair(3) * fpair(2)
1840 tau(9) = tau(9) - dpair(3) * fpair(3)
1841
1842 end subroutine add_stress_tensor
1843
1844 end module doForces