ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 3397
Committed: Tue May 27 16:39:06 2008 UTC (16 years, 1 month ago) by chuckv
File size: 59240 byte(s)
Log Message:
Checking in changes for Hefland moment calculations

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