ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 3129
Committed: Fri Apr 20 18:15:48 2007 UTC (17 years, 3 months ago) by chrisfen
File size: 57821 byte(s)
Log Message:
SF Lennard-Jones was added for everyones' enjoyment.  The behavior is tethered to the electrostaticSummationMethod keyword.

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.87 2007-04-20 18:15:46 chrisfen Exp $, $Date: 2007-04-20 18:15:46 $, $Name: not supported by cvs2svn $, $Revision: 1.87 $
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 suttonchen
66 use status
67 #ifdef IS_MPI
68 use mpiSimulation
69 #endif
70
71 implicit none
72 PRIVATE
73
74 #define __FORTRAN90
75 #include "UseTheForce/fCutoffPolicy.h"
76 #include "UseTheForce/DarkSide/fInteractionMap.h"
77 #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
78
79 INTEGER, PARAMETER:: PREPAIR_LOOP = 1
80 INTEGER, PARAMETER:: PAIR_LOOP = 2
81
82 logical, save :: haveNeighborList = .false.
83 logical, save :: haveSIMvariables = .false.
84 logical, save :: haveSaneForceField = .false.
85 logical, save :: haveInteractionHash = .false.
86 logical, save :: haveGtypeCutoffMap = .false.
87 logical, save :: haveDefaultCutoffs = .false.
88 logical, save :: haveSkinThickness = .false.
89 logical, save :: haveElectrostaticSummationMethod = .false.
90 logical, save :: haveCutoffPolicy = .false.
91 logical, save :: VisitCutoffsAfterComputing = .false.
92 logical, save :: do_box_dipole = .false.
93
94 logical, save :: FF_uses_DirectionalAtoms
95 logical, save :: FF_uses_Dipoles
96 logical, save :: FF_uses_GayBerne
97 logical, save :: FF_uses_EAM
98 logical, save :: FF_uses_SC
99 logical, save :: FF_uses_MEAM
100
101
102 logical, save :: SIM_uses_DirectionalAtoms
103 logical, save :: SIM_uses_EAM
104 logical, save :: SIM_uses_SC
105 logical, save :: SIM_uses_MEAM
106 logical, save :: SIM_requires_postpair_calc
107 logical, save :: SIM_requires_prepair_calc
108 logical, save :: SIM_uses_PBC
109 logical, save :: SIM_uses_AtomicVirial
110
111 integer, save :: electrostaticSummationMethod
112 integer, save :: cutoffPolicy = TRADITIONAL_CUTOFF_POLICY
113
114 real(kind=dp), save :: defaultRcut, defaultRsw, largestRcut
115 real(kind=dp), save :: skinThickness
116 logical, save :: defaultDoShiftPot
117 logical, save :: defaultDoShiftFrc
118
119 public :: init_FF
120 public :: setCutoffs
121 public :: cWasLame
122 public :: setElectrostaticMethod
123 public :: setBoxDipole
124 public :: getBoxDipole
125 public :: setCutoffPolicy
126 public :: setSkinThickness
127 public :: do_force_loop
128
129 #ifdef PROFILE
130 public :: getforcetime
131 real, save :: forceTime = 0
132 real :: forceTimeInitial, forceTimeFinal
133 integer :: nLoops
134 #endif
135
136 !! Variables for cutoff mapping and interaction mapping
137 ! Bit hash to determine pair-pair interactions.
138 integer, dimension(:,:), allocatable :: InteractionHash
139 real(kind=dp), dimension(:), allocatable :: atypeMaxCutoff
140 real(kind=dp), dimension(:), allocatable, target :: groupMaxCutoffRow
141 real(kind=dp), dimension(:), pointer :: groupMaxCutoffCol
142
143 integer, dimension(:), allocatable, target :: groupToGtypeRow
144 integer, dimension(:), pointer :: groupToGtypeCol => null()
145
146 real(kind=dp), dimension(:), allocatable,target :: gtypeMaxCutoffRow
147 real(kind=dp), dimension(:), pointer :: gtypeMaxCutoffCol
148 type ::gtypeCutoffs
149 real(kind=dp) :: rcut
150 real(kind=dp) :: rcutsq
151 real(kind=dp) :: rlistsq
152 end type gtypeCutoffs
153 type(gtypeCutoffs), dimension(:,:), allocatable :: gtypeCutoffMap
154
155 real(kind=dp), dimension(3) :: boxDipole
156
157 contains
158
159 subroutine createInteractionHash()
160 integer :: nAtypes
161 integer :: i
162 integer :: j
163 integer :: iHash
164 !! Test Types
165 logical :: i_is_LJ
166 logical :: i_is_Elect
167 logical :: i_is_Sticky
168 logical :: i_is_StickyP
169 logical :: i_is_GB
170 logical :: i_is_EAM
171 logical :: i_is_Shape
172 logical :: i_is_SC
173 logical :: i_is_MEAM
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 logical :: j_is_MEAM
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 call getElementProperty(atypes, i, "is_MEAM", i_is_MEAM)
221
222 do j = i, nAtypes
223
224 iHash = 0
225 myRcut = 0.0_dp
226
227 call getElementProperty(atypes, j, "is_LennardJones", j_is_LJ)
228 call getElementProperty(atypes, j, "is_Electrostatic", j_is_Elect)
229 call getElementProperty(atypes, j, "is_Sticky", j_is_Sticky)
230 call getElementProperty(atypes, j, "is_StickyPower", j_is_StickyP)
231 call getElementProperty(atypes, j, "is_GayBerne", j_is_GB)
232 call getElementProperty(atypes, j, "is_EAM", j_is_EAM)
233 call getElementProperty(atypes, j, "is_Shape", j_is_Shape)
234 call getElementProperty(atypes, j, "is_SC", j_is_SC)
235 call getElementProperty(atypes, j, "is_MEAM", j_is_MEAM)
236
237 if (i_is_LJ .and. j_is_LJ) then
238 iHash = ior(iHash, LJ_PAIR)
239 endif
240
241 if (i_is_Elect .and. j_is_Elect) then
242 iHash = ior(iHash, ELECTROSTATIC_PAIR)
243 endif
244
245 if (i_is_Sticky .and. j_is_Sticky) then
246 iHash = ior(iHash, STICKY_PAIR)
247 endif
248
249 if (i_is_StickyP .and. j_is_StickyP) then
250 iHash = ior(iHash, STICKYPOWER_PAIR)
251 endif
252
253 if (i_is_EAM .and. j_is_EAM) then
254 iHash = ior(iHash, EAM_PAIR)
255 endif
256
257 if (i_is_SC .and. j_is_SC) then
258 iHash = ior(iHash, SC_PAIR)
259 endif
260
261 if (i_is_GB .and. j_is_GB) iHash = ior(iHash, GAYBERNE_PAIR)
262 if (i_is_GB .and. j_is_LJ) iHash = ior(iHash, GAYBERNE_LJ)
263 if (i_is_LJ .and. j_is_GB) iHash = ior(iHash, GAYBERNE_LJ)
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 set_switch(defaultRsw, defaultRcut)
607 call setHmatDangerousRcutValue(defaultRcut)
608
609 haveDefaultCutoffs = .true.
610 haveGtypeCutoffMap = .false.
611
612 end subroutine setCutoffs
613
614 subroutine cWasLame()
615
616 VisitCutoffsAfterComputing = .true.
617 return
618
619 end subroutine cWasLame
620
621 subroutine setCutoffPolicy(cutPolicy)
622
623 integer, intent(in) :: cutPolicy
624
625 cutoffPolicy = cutPolicy
626 haveCutoffPolicy = .true.
627 haveGtypeCutoffMap = .false.
628
629 end subroutine setCutoffPolicy
630
631 subroutine setBoxDipole()
632
633 do_box_dipole = .true.
634
635 end subroutine setBoxDipole
636
637 subroutine getBoxDipole( box_dipole )
638
639 real(kind=dp), intent(inout), dimension(3) :: box_dipole
640
641 box_dipole = boxDipole
642
643 end subroutine getBoxDipole
644
645 subroutine setElectrostaticMethod( thisESM )
646
647 integer, intent(in) :: thisESM
648
649 electrostaticSummationMethod = thisESM
650 haveElectrostaticSummationMethod = .true.
651
652 end subroutine setElectrostaticMethod
653
654 subroutine setSkinThickness( thisSkin )
655
656 real(kind=dp), intent(in) :: thisSkin
657
658 skinThickness = thisSkin
659 haveSkinThickness = .true.
660 haveGtypeCutoffMap = .false.
661
662 end subroutine setSkinThickness
663
664 subroutine setSimVariables()
665 SIM_uses_DirectionalAtoms = SimUsesDirectionalAtoms()
666 SIM_uses_EAM = SimUsesEAM()
667 SIM_requires_postpair_calc = SimRequiresPostpairCalc()
668 SIM_requires_prepair_calc = SimRequiresPrepairCalc()
669 SIM_uses_PBC = SimUsesPBC()
670 SIM_uses_SC = SimUsesSC()
671 SIM_uses_AtomicVirial = SimUsesAtomicVirial()
672
673 haveSIMvariables = .true.
674
675 return
676 end subroutine setSimVariables
677
678 subroutine doReadyCheck(error)
679 integer, intent(out) :: error
680 integer :: myStatus
681
682 error = 0
683
684 if (.not. haveInteractionHash) then
685 call createInteractionHash()
686 endif
687
688 if (.not. haveGtypeCutoffMap) then
689 call createGtypeCutoffMap()
690 endif
691
692 if (VisitCutoffsAfterComputing) then
693 call set_switch(largestRcut, largestRcut)
694 call setHmatDangerousRcutValue(largestRcut)
695 call setCutoffEAM(largestRcut)
696 call setCutoffSC(largestRcut)
697 VisitCutoffsAfterComputing = .false.
698 endif
699
700 if (.not. haveSIMvariables) then
701 call setSimVariables()
702 endif
703
704 if (.not. haveNeighborList) then
705 write(default_error, *) 'neighbor list has not been initialized in doForces!'
706 error = -1
707 return
708 end if
709
710 if (.not. haveSaneForceField) then
711 write(default_error, *) 'Force Field is not sane in doForces!'
712 error = -1
713 return
714 end if
715
716 #ifdef IS_MPI
717 if (.not. isMPISimSet()) then
718 write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
719 error = -1
720 return
721 endif
722 #endif
723 return
724 end subroutine doReadyCheck
725
726
727 subroutine init_FF(thisStat)
728
729 integer, intent(out) :: thisStat
730 integer :: my_status, nMatches
731 integer, pointer :: MatchList(:) => null()
732
733 !! assume things are copacetic, unless they aren't
734 thisStat = 0
735
736 !! init_FF is called *after* all of the atom types have been
737 !! defined in atype_module using the new_atype subroutine.
738 !!
739 !! this will scan through the known atypes and figure out what
740 !! interactions are used by the force field.
741
742 FF_uses_DirectionalAtoms = .false.
743 FF_uses_Dipoles = .false.
744 FF_uses_GayBerne = .false.
745 FF_uses_EAM = .false.
746 FF_uses_SC = .false.
747
748 call getMatchingElementList(atypes, "is_Directional", .true., &
749 nMatches, MatchList)
750 if (nMatches .gt. 0) FF_uses_DirectionalAtoms = .true.
751
752 call getMatchingElementList(atypes, "is_Dipole", .true., &
753 nMatches, MatchList)
754 if (nMatches .gt. 0) FF_uses_Dipoles = .true.
755
756 call getMatchingElementList(atypes, "is_GayBerne", .true., &
757 nMatches, MatchList)
758 if (nMatches .gt. 0) FF_uses_GayBerne = .true.
759
760 call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList)
761 if (nMatches .gt. 0) FF_uses_EAM = .true.
762
763 call getMatchingElementList(atypes, "is_SC", .true., nMatches, MatchList)
764 if (nMatches .gt. 0) FF_uses_SC = .true.
765
766
767 haveSaneForceField = .true.
768
769 if (FF_uses_EAM) then
770 call init_EAM_FF(my_status)
771 if (my_status /= 0) then
772 write(default_error, *) "init_EAM_FF returned a bad status"
773 thisStat = -1
774 haveSaneForceField = .false.
775 return
776 end if
777 endif
778
779 if (.not. haveNeighborList) then
780 !! Create neighbor lists
781 call expandNeighborList(nLocal, my_status)
782 if (my_Status /= 0) then
783 write(default_error,*) "SimSetup: ExpandNeighborList returned error."
784 thisStat = -1
785 return
786 endif
787 haveNeighborList = .true.
788 endif
789
790 end subroutine init_FF
791
792
793 !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
794 !------------------------------------------------------------->
795 subroutine do_force_loop(q, q_group, A, eFrame, f, t, tau, pot, &
796 do_pot_c, do_stress_c, error)
797 !! Position array provided by C, dimensioned by getNlocal
798 real ( kind = dp ), dimension(3, nLocal) :: q
799 !! molecular center-of-mass position array
800 real ( kind = dp ), dimension(3, nGroups) :: q_group
801 !! Rotation Matrix for each long range particle in simulation.
802 real( kind = dp), dimension(9, nLocal) :: A
803 !! Unit vectors for dipoles (lab frame)
804 real( kind = dp ), dimension(9,nLocal) :: eFrame
805 !! Force array provided by C, dimensioned by getNlocal
806 real ( kind = dp ), dimension(3,nLocal) :: f
807 !! Torsion array provided by C, dimensioned by getNlocal
808 real( kind = dp ), dimension(3,nLocal) :: t
809
810 !! Stress Tensor
811 real( kind = dp), dimension(9) :: tau
812 real ( kind = dp ),dimension(LR_POT_TYPES) :: pot
813 logical ( kind = 2) :: do_pot_c, do_stress_c
814 logical :: do_pot
815 logical :: do_stress
816 logical :: in_switching_region
817 #ifdef IS_MPI
818 real( kind = DP ), dimension(LR_POT_TYPES) :: pot_local
819 integer :: nAtomsInRow
820 integer :: nAtomsInCol
821 integer :: nprocs
822 integer :: nGroupsInRow
823 integer :: nGroupsInCol
824 #endif
825 integer :: natoms
826 logical :: update_nlist
827 integer :: i, j, jstart, jend, jnab
828 integer :: istart, iend
829 integer :: ia, jb, atom1, atom2
830 integer :: nlist
831 real( kind = DP ) :: ratmsq, rgrpsq, rgrp, rag, vpair, vij
832 real( kind = DP ) :: sw, dswdr, swderiv, mf
833 real( kind = DP ) :: rVal
834 real(kind=dp),dimension(3) :: d_atm, d_grp, fpair, fij, fg, dag
835 real(kind=dp) :: rfpot, mu_i
836 real(kind=dp):: rCut
837 integer :: me_i, me_j, n_in_i, n_in_j
838 logical :: is_dp_i
839 integer :: neighborListSize
840 integer :: listerror, error
841 integer :: localError
842 integer :: propPack_i, propPack_j
843 integer :: loopStart, loopEnd, loop
844 integer :: iHash
845 integer :: i1
846
847 !! the variables for the box dipole moment
848 #ifdef IS_MPI
849 integer :: pChgCount_local
850 integer :: nChgCount_local
851 real(kind=dp) :: pChg_local
852 real(kind=dp) :: nChg_local
853 real(kind=dp), dimension(3) :: pChgPos_local
854 real(kind=dp), dimension(3) :: nChgPos_local
855 real(kind=dp), dimension(3) :: dipVec_local
856 #endif
857 integer :: pChgCount
858 integer :: nChgCount
859 real(kind=dp) :: pChg
860 real(kind=dp) :: nChg
861 real(kind=dp) :: chg_value
862 real(kind=dp), dimension(3) :: pChgPos
863 real(kind=dp), dimension(3) :: nChgPos
864 real(kind=dp), dimension(3) :: dipVec
865 real(kind=dp), dimension(3) :: chgVec
866
867 !! initialize box dipole variables
868 if (do_box_dipole) then
869 #ifdef IS_MPI
870 pChg_local = 0.0_dp
871 nChg_local = 0.0_dp
872 pChgCount_local = 0
873 nChgCount_local = 0
874 do i=1, 3
875 pChgPos_local = 0.0_dp
876 nChgPos_local = 0.0_dp
877 dipVec_local = 0.0_dp
878 enddo
879 #endif
880 pChg = 0.0_dp
881 nChg = 0.0_dp
882 pChgCount = 0
883 nChgCount = 0
884 chg_value = 0.0_dp
885
886 do i=1, 3
887 pChgPos(i) = 0.0_dp
888 nChgPos(i) = 0.0_dp
889 dipVec(i) = 0.0_dp
890 chgVec(i) = 0.0_dp
891 boxDipole(i) = 0.0_dp
892 enddo
893 endif
894
895 !! initialize local variables
896
897 #ifdef IS_MPI
898 pot_local = 0.0_dp
899 nAtomsInRow = getNatomsInRow(plan_atom_row)
900 nAtomsInCol = getNatomsInCol(plan_atom_col)
901 nGroupsInRow = getNgroupsInRow(plan_group_row)
902 nGroupsInCol = getNgroupsInCol(plan_group_col)
903 #else
904 natoms = nlocal
905 #endif
906
907 call doReadyCheck(localError)
908 if ( localError .ne. 0 ) then
909 call handleError("do_force_loop", "Not Initialized")
910 error = -1
911 return
912 end if
913 call zero_work_arrays()
914
915 do_pot = do_pot_c
916 do_stress = do_stress_c
917
918 ! Gather all information needed by all force loops:
919
920 #ifdef IS_MPI
921
922 call gather(q, q_Row, plan_atom_row_3d)
923 call gather(q, q_Col, plan_atom_col_3d)
924
925 call gather(q_group, q_group_Row, plan_group_row_3d)
926 call gather(q_group, q_group_Col, plan_group_col_3d)
927
928 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
929 call gather(eFrame, eFrame_Row, plan_atom_row_rotation)
930 call gather(eFrame, eFrame_Col, plan_atom_col_rotation)
931
932 call gather(A, A_Row, plan_atom_row_rotation)
933 call gather(A, A_Col, plan_atom_col_rotation)
934 endif
935
936 #endif
937
938 !! Begin force loop timing:
939 #ifdef PROFILE
940 call cpu_time(forceTimeInitial)
941 nloops = nloops + 1
942 #endif
943
944 loopEnd = PAIR_LOOP
945 if (FF_RequiresPrepairCalc() .and. SIM_requires_prepair_calc) then
946 loopStart = PREPAIR_LOOP
947 else
948 loopStart = PAIR_LOOP
949 endif
950
951 do loop = loopStart, loopEnd
952
953 ! See if we need to update neighbor lists
954 ! (but only on the first time through):
955 if (loop .eq. loopStart) then
956 #ifdef IS_MPI
957 call checkNeighborList(nGroupsInRow, q_group_row, skinThickness, &
958 update_nlist)
959 #else
960 call checkNeighborList(nGroups, q_group, skinThickness, &
961 update_nlist)
962 #endif
963 endif
964
965 if (update_nlist) then
966 !! save current configuration and construct neighbor list
967 #ifdef IS_MPI
968 call saveNeighborList(nGroupsInRow, q_group_row)
969 #else
970 call saveNeighborList(nGroups, q_group)
971 #endif
972 neighborListSize = size(list)
973 nlist = 0
974 endif
975
976 istart = 1
977 #ifdef IS_MPI
978 iend = nGroupsInRow
979 #else
980 iend = nGroups - 1
981 #endif
982 outer: do i = istart, iend
983
984 if (update_nlist) point(i) = nlist + 1
985
986 n_in_i = groupStartRow(i+1) - groupStartRow(i)
987
988 if (update_nlist) then
989 #ifdef IS_MPI
990 jstart = 1
991 jend = nGroupsInCol
992 #else
993 jstart = i+1
994 jend = nGroups
995 #endif
996 else
997 jstart = point(i)
998 jend = point(i+1) - 1
999 ! make sure group i has neighbors
1000 if (jstart .gt. jend) cycle outer
1001 endif
1002
1003 do jnab = jstart, jend
1004 if (update_nlist) then
1005 j = jnab
1006 else
1007 j = list(jnab)
1008 endif
1009
1010 #ifdef IS_MPI
1011 me_j = atid_col(j)
1012 call get_interatomic_vector(q_group_Row(:,i), &
1013 q_group_Col(:,j), d_grp, rgrpsq)
1014 #else
1015 me_j = atid(j)
1016 call get_interatomic_vector(q_group(:,i), &
1017 q_group(:,j), d_grp, rgrpsq)
1018 #endif
1019
1020 if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rListsq) then
1021 if (update_nlist) then
1022 nlist = nlist + 1
1023
1024 if (nlist > neighborListSize) then
1025 #ifdef IS_MPI
1026 call expandNeighborList(nGroupsInRow, listerror)
1027 #else
1028 call expandNeighborList(nGroups, listerror)
1029 #endif
1030 if (listerror /= 0) then
1031 error = -1
1032 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
1033 return
1034 end if
1035 neighborListSize = size(list)
1036 endif
1037
1038 list(nlist) = j
1039 endif
1040
1041 if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCutsq) then
1042
1043 rCut = gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCut
1044 if (loop .eq. PAIR_LOOP) then
1045 vij = 0.0_dp
1046 fij(1) = 0.0_dp
1047 fij(2) = 0.0_dp
1048 fij(3) = 0.0_dp
1049 endif
1050
1051 call get_switch(rgrpsq, sw, dswdr,rgrp, in_switching_region)
1052
1053 n_in_j = groupStartCol(j+1) - groupStartCol(j)
1054
1055 do ia = groupStartRow(i), groupStartRow(i+1)-1
1056
1057 atom1 = groupListRow(ia)
1058
1059 inner: do jb = groupStartCol(j), groupStartCol(j+1)-1
1060
1061 atom2 = groupListCol(jb)
1062
1063 if (skipThisPair(atom1, atom2)) cycle inner
1064
1065 if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
1066 d_atm(1) = d_grp(1)
1067 d_atm(2) = d_grp(2)
1068 d_atm(3) = d_grp(3)
1069 ratmsq = rgrpsq
1070 else
1071 #ifdef IS_MPI
1072 call get_interatomic_vector(q_Row(:,atom1), &
1073 q_Col(:,atom2), d_atm, ratmsq)
1074 #else
1075 call get_interatomic_vector(q(:,atom1), &
1076 q(:,atom2), d_atm, ratmsq)
1077 #endif
1078 endif
1079
1080 if (loop .eq. PREPAIR_LOOP) then
1081 #ifdef IS_MPI
1082 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
1083 rgrpsq, d_grp, rCut, do_pot, do_stress, &
1084 eFrame, A, f, t, pot_local)
1085 #else
1086 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
1087 rgrpsq, d_grp, rCut, do_pot, do_stress, &
1088 eFrame, A, f, t, pot)
1089 #endif
1090 else
1091 #ifdef IS_MPI
1092 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
1093 do_pot, eFrame, A, f, t, pot_local, vpair, &
1094 fpair, d_grp, rgrp, rCut)
1095 #else
1096 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
1097 do_pot, eFrame, A, f, t, pot, vpair, fpair, &
1098 d_grp, rgrp, rCut)
1099 #endif
1100 vij = vij + vpair
1101 fij(1) = fij(1) + fpair(1)
1102 fij(2) = fij(2) + fpair(2)
1103 fij(3) = fij(3) + fpair(3)
1104 if (do_stress) then
1105 call add_stress_tensor(d_atm, fpair, tau)
1106 endif
1107 endif
1108 enddo inner
1109 enddo
1110
1111 if (loop .eq. PAIR_LOOP) then
1112 if (in_switching_region) then
1113 swderiv = vij*dswdr/rgrp
1114 fij(1) = fij(1) + swderiv*d_grp(1)
1115 fij(2) = fij(2) + swderiv*d_grp(2)
1116 fij(3) = fij(3) + swderiv*d_grp(3)
1117
1118 do ia=groupStartRow(i), groupStartRow(i+1)-1
1119 atom1=groupListRow(ia)
1120 mf = mfactRow(atom1)
1121 ! fg is the force on atom ia due to cutoff group's
1122 ! presence in switching region
1123 fg = swderiv*d_grp*mf
1124 #ifdef IS_MPI
1125 f_Row(1,atom1) = f_Row(1,atom1) + fg(1)
1126 f_Row(2,atom1) = f_Row(2,atom1) + fg(2)
1127 f_Row(3,atom1) = f_Row(3,atom1) + fg(3)
1128 #else
1129 f(1,atom1) = f(1,atom1) + fg(1)
1130 f(2,atom1) = f(2,atom1) + fg(2)
1131 f(3,atom1) = f(3,atom1) + fg(3)
1132 #endif
1133 if (n_in_i .gt. 1) then
1134 if (do_stress.and.SIM_uses_AtomicVirial) then
1135 ! find the distance between the atom and the center of
1136 ! the cutoff group:
1137 #ifdef IS_MPI
1138 call get_interatomic_vector(q_Row(:,atom1), &
1139 q_group_Row(:,i), dag, rag)
1140 #else
1141 call get_interatomic_vector(q(:,atom1), &
1142 q_group(:,i), dag, rag)
1143 #endif
1144 call add_stress_tensor(dag,fg,tau)
1145 endif
1146 endif
1147 enddo
1148
1149 do jb=groupStartCol(j), groupStartCol(j+1)-1
1150 atom2=groupListCol(jb)
1151 mf = mfactCol(atom2)
1152 ! fg is the force on atom jb due to cutoff group's
1153 ! presence in switching region
1154 fg = -swderiv*d_grp*mf
1155 #ifdef IS_MPI
1156 f_Col(1,atom2) = f_Col(1,atom2) + fg(1)
1157 f_Col(2,atom2) = f_Col(2,atom2) + fg(2)
1158 f_Col(3,atom2) = f_Col(3,atom2) + fg(3)
1159 #else
1160 f(1,atom2) = f(1,atom2) + fg(1)
1161 f(2,atom2) = f(2,atom2) + fg(2)
1162 f(3,atom2) = f(3,atom2) + fg(3)
1163 #endif
1164 if (n_in_j .gt. 1) then
1165 if (do_stress.and.SIM_uses_AtomicVirial) then
1166 ! find the distance between the atom and the center of
1167 ! the cutoff group:
1168 #ifdef IS_MPI
1169 call get_interatomic_vector(q_Col(:,atom2), &
1170 q_group_Col(:,j), dag, rag)
1171 #else
1172 call get_interatomic_vector(q(:,atom2), &
1173 q_group(:,j), dag, rag)
1174 #endif
1175 call add_stress_tensor(dag,fg,tau)
1176 endif
1177 endif
1178 enddo
1179 endif
1180 endif
1181 endif
1182 endif
1183 enddo
1184
1185 enddo outer
1186
1187 if (update_nlist) then
1188 #ifdef IS_MPI
1189 point(nGroupsInRow + 1) = nlist + 1
1190 #else
1191 point(nGroups) = nlist + 1
1192 #endif
1193 if (loop .eq. PREPAIR_LOOP) then
1194 ! we just did the neighbor list update on the first
1195 ! pass, so we don't need to do it
1196 ! again on the second pass
1197 update_nlist = .false.
1198 endif
1199 endif
1200
1201 if (loop .eq. PREPAIR_LOOP) then
1202 call do_preforce(nlocal, pot)
1203 endif
1204
1205 enddo
1206
1207 !! Do timing
1208 #ifdef PROFILE
1209 call cpu_time(forceTimeFinal)
1210 forceTime = forceTime + forceTimeFinal - forceTimeInitial
1211 #endif
1212
1213 #ifdef IS_MPI
1214 !!distribute forces
1215
1216 f_temp = 0.0_dp
1217 call scatter(f_Row,f_temp,plan_atom_row_3d)
1218 do i = 1,nlocal
1219 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1220 end do
1221
1222 f_temp = 0.0_dp
1223 call scatter(f_Col,f_temp,plan_atom_col_3d)
1224 do i = 1,nlocal
1225 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1226 end do
1227
1228 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
1229 t_temp = 0.0_dp
1230 call scatter(t_Row,t_temp,plan_atom_row_3d)
1231 do i = 1,nlocal
1232 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1233 end do
1234 t_temp = 0.0_dp
1235 call scatter(t_Col,t_temp,plan_atom_col_3d)
1236
1237 do i = 1,nlocal
1238 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1239 end do
1240 endif
1241
1242 if (do_pot) then
1243 ! scatter/gather pot_row into the members of my column
1244 do i = 1,LR_POT_TYPES
1245 call scatter(pot_Row(i,:), pot_Temp(i,:), plan_atom_row)
1246 end do
1247 ! scatter/gather pot_local into all other procs
1248 ! add resultant to get total pot
1249 do i = 1, nlocal
1250 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES) &
1251 + pot_Temp(1:LR_POT_TYPES,i)
1252 enddo
1253
1254 pot_Temp = 0.0_DP
1255 do i = 1,LR_POT_TYPES
1256 call scatter(pot_Col(i,:), pot_Temp(i,:), plan_atom_col)
1257 end do
1258 do i = 1, nlocal
1259 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES)&
1260 + pot_Temp(1:LR_POT_TYPES,i)
1261 enddo
1262
1263 endif
1264 #endif
1265
1266 if (SIM_requires_postpair_calc) then
1267 do i = 1, nlocal
1268
1269 ! we loop only over the local atoms, so we don't need row and column
1270 ! lookups for the types
1271
1272 me_i = atid(i)
1273
1274 ! is the atom electrostatic? See if it would have an
1275 ! electrostatic interaction with itself
1276 iHash = InteractionHash(me_i,me_i)
1277
1278 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1279 #ifdef IS_MPI
1280 call self_self(i, eFrame, pot_local(ELECTROSTATIC_POT), &
1281 t, do_pot)
1282 #else
1283 call self_self(i, eFrame, pot(ELECTROSTATIC_POT), &
1284 t, do_pot)
1285 #endif
1286 endif
1287
1288
1289 if (electrostaticSummationMethod.eq.REACTION_FIELD) then
1290
1291 ! loop over the excludes to accumulate RF stuff we've
1292 ! left out of the normal pair loop
1293
1294 do i1 = 1, nSkipsForAtom(i)
1295 j = skipsForAtom(i, i1)
1296
1297 ! prevent overcounting of the skips
1298 if (i.lt.j) then
1299 call get_interatomic_vector(q(:,i), q(:,j), d_atm, ratmsq)
1300 rVal = sqrt(ratmsq)
1301 call get_switch(ratmsq, sw, dswdr, rVal,in_switching_region)
1302 #ifdef IS_MPI
1303 call rf_self_excludes(i, j, sw, eFrame, d_atm, rVal, &
1304 vpair, pot_local(ELECTROSTATIC_POT), f, t, do_pot)
1305 #else
1306 call rf_self_excludes(i, j, sw, eFrame, d_atm, rVal, &
1307 vpair, pot(ELECTROSTATIC_POT), f, t, do_pot)
1308 #endif
1309 endif
1310 enddo
1311 endif
1312
1313 if (do_box_dipole) then
1314 #ifdef IS_MPI
1315 call accumulate_box_dipole(i, eFrame, q(:,i), pChg_local, &
1316 nChg_local, pChgPos_local, nChgPos_local, dipVec_local, &
1317 pChgCount_local, nChgCount_local)
1318 #else
1319 call accumulate_box_dipole(i, eFrame, q(:,i), pChg, nChg, &
1320 pChgPos, nChgPos, dipVec, pChgCount, nChgCount)
1321 #endif
1322 endif
1323 enddo
1324 endif
1325
1326 #ifdef IS_MPI
1327 if (do_pot) then
1328 #ifdef SINGLE_PRECISION
1329 call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_real,mpi_sum, &
1330 mpi_comm_world,mpi_err)
1331 #else
1332 call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_double_precision, &
1333 mpi_sum, mpi_comm_world,mpi_err)
1334 #endif
1335 endif
1336
1337 if (do_box_dipole) then
1338
1339 #ifdef SINGLE_PRECISION
1340 call mpi_allreduce(pChg_local, pChg, 1, mpi_real, mpi_sum, &
1341 mpi_comm_world, mpi_err)
1342 call mpi_allreduce(nChg_local, nChg, 1, mpi_real, mpi_sum, &
1343 mpi_comm_world, mpi_err)
1344 call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer, mpi_sum,&
1345 mpi_comm_world, mpi_err)
1346 call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer, mpi_sum,&
1347 mpi_comm_world, mpi_err)
1348 call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_real, mpi_sum, &
1349 mpi_comm_world, mpi_err)
1350 call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_real, mpi_sum, &
1351 mpi_comm_world, mpi_err)
1352 call mpi_allreduce(dipVec_local, dipVec, 3, mpi_real, mpi_sum, &
1353 mpi_comm_world, mpi_err)
1354 #else
1355 call mpi_allreduce(pChg_local, pChg, 1, mpi_double_precision, mpi_sum, &
1356 mpi_comm_world, mpi_err)
1357 call mpi_allreduce(nChg_local, nChg, 1, mpi_double_precision, mpi_sum, &
1358 mpi_comm_world, mpi_err)
1359 call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer,&
1360 mpi_sum, mpi_comm_world, mpi_err)
1361 call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer,&
1362 mpi_sum, mpi_comm_world, mpi_err)
1363 call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_double_precision, &
1364 mpi_sum, mpi_comm_world, mpi_err)
1365 call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_double_precision, &
1366 mpi_sum, mpi_comm_world, mpi_err)
1367 call mpi_allreduce(dipVec_local, dipVec, 3, mpi_double_precision, &
1368 mpi_sum, mpi_comm_world, mpi_err)
1369 #endif
1370
1371 endif
1372
1373 #endif
1374
1375 if (do_box_dipole) then
1376 ! first load the accumulated dipole moment (if dipoles were present)
1377 boxDipole(1) = dipVec(1)
1378 boxDipole(2) = dipVec(2)
1379 boxDipole(3) = dipVec(3)
1380
1381 ! now include the dipole moment due to charges
1382 ! use the lesser of the positive and negative charge totals
1383 if (nChg .le. pChg) then
1384 chg_value = nChg
1385 else
1386 chg_value = pChg
1387 endif
1388
1389 ! find the average positions
1390 if (pChgCount .gt. 0 .and. nChgCount .gt. 0) then
1391 pChgPos = pChgPos / pChgCount
1392 nChgPos = nChgPos / nChgCount
1393 endif
1394
1395 ! dipole is from the negative to the positive (physics notation)
1396 chgVec(1) = pChgPos(1) - nChgPos(1)
1397 chgVec(2) = pChgPos(2) - nChgPos(2)
1398 chgVec(3) = pChgPos(3) - nChgPos(3)
1399
1400 boxDipole(1) = boxDipole(1) + chgVec(1) * chg_value
1401 boxDipole(2) = boxDipole(2) + chgVec(2) * chg_value
1402 boxDipole(3) = boxDipole(3) + chgVec(3) * chg_value
1403
1404 endif
1405
1406 end subroutine do_force_loop
1407
1408 subroutine do_pair(i, j, rijsq, d, sw, do_pot, &
1409 eFrame, A, f, t, pot, vpair, fpair, d_grp, r_grp, rCut)
1410
1411 real( kind = dp ) :: vpair, sw
1412 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1413 real( kind = dp ), dimension(3) :: fpair
1414 real( kind = dp ), dimension(nLocal) :: mfact
1415 real( kind = dp ), dimension(9,nLocal) :: eFrame
1416 real( kind = dp ), dimension(9,nLocal) :: A
1417 real( kind = dp ), dimension(3,nLocal) :: f
1418 real( kind = dp ), dimension(3,nLocal) :: t
1419
1420 logical, intent(inout) :: do_pot
1421 integer, intent(in) :: i, j
1422 real ( kind = dp ), intent(inout) :: rijsq
1423 real ( kind = dp ), intent(inout) :: r_grp
1424 real ( kind = dp ), intent(inout) :: d(3)
1425 real ( kind = dp ), intent(inout) :: d_grp(3)
1426 real ( kind = dp ), intent(inout) :: rCut
1427 real ( kind = dp ) :: r
1428 real ( kind = dp ) :: a_k, b_k, c_k, d_k, dx
1429 integer :: me_i, me_j
1430 integer :: k
1431
1432 integer :: iHash
1433
1434 r = sqrt(rijsq)
1435
1436 vpair = 0.0_dp
1437 fpair(1:3) = 0.0_dp
1438
1439 #ifdef IS_MPI
1440 me_i = atid_row(i)
1441 me_j = atid_col(j)
1442 #else
1443 me_i = atid(i)
1444 me_j = atid(j)
1445 #endif
1446
1447 iHash = InteractionHash(me_i, me_j)
1448
1449 if ( iand(iHash, LJ_PAIR).ne.0 ) then
1450 call do_lj_pair(i, j, d, r, rijsq, rcut, sw, vpair, fpair, &
1451 pot(VDW_POT), f, do_pot)
1452 endif
1453
1454 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1455 call doElectrostaticPair(i, j, d, r, rijsq, rcut, sw, vpair, fpair, &
1456 pot(ELECTROSTATIC_POT), eFrame, f, t, do_pot)
1457 endif
1458
1459 if ( iand(iHash, STICKY_PAIR).ne.0 ) then
1460 call do_sticky_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1461 pot(HB_POT), A, f, t, do_pot)
1462 endif
1463
1464 if ( iand(iHash, STICKYPOWER_PAIR).ne.0 ) then
1465 call do_sticky_power_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1466 pot(HB_POT), A, f, t, do_pot)
1467 endif
1468
1469 if ( iand(iHash, GAYBERNE_PAIR).ne.0 ) then
1470 call do_gb_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1471 pot(VDW_POT), A, f, t, do_pot)
1472 endif
1473
1474 if ( iand(iHash, GAYBERNE_LJ).ne.0 ) then
1475 call do_gb_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1476 pot(VDW_POT), A, f, t, do_pot)
1477 endif
1478
1479 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1480 call do_eam_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1481 pot(METALLIC_POT), f, do_pot)
1482 endif
1483
1484 if ( iand(iHash, SHAPE_PAIR).ne.0 ) then
1485 call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1486 pot(VDW_POT), A, f, t, do_pot)
1487 endif
1488
1489 if ( iand(iHash, SHAPE_LJ).ne.0 ) then
1490 call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1491 pot(VDW_POT), A, f, t, do_pot)
1492 endif
1493
1494 if ( iand(iHash, SC_PAIR).ne.0 ) then
1495 call do_SC_pair(i, j, d, r, rijsq, rcut, sw, vpair, fpair, &
1496 pot(METALLIC_POT), f, do_pot)
1497 endif
1498
1499 end subroutine do_pair
1500
1501 subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, rCut, &
1502 do_pot, do_stress, eFrame, A, f, t, pot)
1503
1504 real( kind = dp ) :: sw
1505 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1506 real( kind = dp ), dimension(9,nLocal) :: eFrame
1507 real (kind=dp), dimension(9,nLocal) :: A
1508 real (kind=dp), dimension(3,nLocal) :: f
1509 real (kind=dp), dimension(3,nLocal) :: t
1510
1511 logical, intent(inout) :: do_pot, do_stress
1512 integer, intent(in) :: i, j
1513 real ( kind = dp ), intent(inout) :: rijsq, rcijsq, rCut
1514 real ( kind = dp ) :: r, rc
1515 real ( kind = dp ), intent(inout) :: d(3), dc(3)
1516
1517 integer :: me_i, me_j, iHash
1518
1519 r = sqrt(rijsq)
1520
1521 #ifdef IS_MPI
1522 me_i = atid_row(i)
1523 me_j = atid_col(j)
1524 #else
1525 me_i = atid(i)
1526 me_j = atid(j)
1527 #endif
1528
1529 iHash = InteractionHash(me_i, me_j)
1530
1531 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1532 call calc_EAM_prepair_rho(i, j, d, r, rijsq)
1533 endif
1534
1535 if ( iand(iHash, SC_PAIR).ne.0 ) then
1536 call calc_SC_prepair_rho(i, j, d, r, rijsq, rcut )
1537 endif
1538
1539 end subroutine do_prepair
1540
1541
1542 subroutine do_preforce(nlocal,pot)
1543 integer :: nlocal
1544 real( kind = dp ),dimension(LR_POT_TYPES) :: pot
1545
1546 if (FF_uses_EAM .and. SIM_uses_EAM) then
1547 call calc_EAM_preforce_Frho(nlocal,pot(METALLIC_POT))
1548 endif
1549 if (FF_uses_SC .and. SIM_uses_SC) then
1550 call calc_SC_preforce_Frho(nlocal,pot(METALLIC_POT))
1551 endif
1552 end subroutine do_preforce
1553
1554
1555 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
1556
1557 real (kind = dp), dimension(3) :: q_i
1558 real (kind = dp), dimension(3) :: q_j
1559 real ( kind = dp ), intent(out) :: r_sq
1560 real( kind = dp ) :: d(3), scaled(3)
1561 integer i
1562
1563 d(1) = q_j(1) - q_i(1)
1564 d(2) = q_j(2) - q_i(2)
1565 d(3) = q_j(3) - q_i(3)
1566
1567 ! Wrap back into periodic box if necessary
1568 if ( SIM_uses_PBC ) then
1569
1570 if( .not.boxIsOrthorhombic ) then
1571 ! calc the scaled coordinates.
1572 ! scaled = matmul(HmatInv, d)
1573
1574 scaled(1) = HmatInv(1,1)*d(1) + HmatInv(1,2)*d(2) + HmatInv(1,3)*d(3)
1575 scaled(2) = HmatInv(2,1)*d(1) + HmatInv(2,2)*d(2) + HmatInv(2,3)*d(3)
1576 scaled(3) = HmatInv(3,1)*d(1) + HmatInv(3,2)*d(2) + HmatInv(3,3)*d(3)
1577
1578 ! wrap the scaled coordinates
1579
1580 scaled(1) = scaled(1) - anint(scaled(1), kind=dp)
1581 scaled(2) = scaled(2) - anint(scaled(2), kind=dp)
1582 scaled(3) = scaled(3) - anint(scaled(3), kind=dp)
1583
1584 ! calc the wrapped real coordinates from the wrapped scaled
1585 ! coordinates
1586 ! d = matmul(Hmat,scaled)
1587 d(1)= Hmat(1,1)*scaled(1) + Hmat(1,2)*scaled(2) + Hmat(1,3)*scaled(3)
1588 d(2)= Hmat(2,1)*scaled(1) + Hmat(2,2)*scaled(2) + Hmat(2,3)*scaled(3)
1589 d(3)= Hmat(3,1)*scaled(1) + Hmat(3,2)*scaled(2) + Hmat(3,3)*scaled(3)
1590
1591 else
1592 ! calc the scaled coordinates.
1593
1594 scaled(1) = d(1) * HmatInv(1,1)
1595 scaled(2) = d(2) * HmatInv(2,2)
1596 scaled(3) = d(3) * HmatInv(3,3)
1597
1598 ! wrap the scaled coordinates
1599
1600 scaled(1) = scaled(1) - anint(scaled(1), kind=dp)
1601 scaled(2) = scaled(2) - anint(scaled(2), kind=dp)
1602 scaled(3) = scaled(3) - anint(scaled(3), kind=dp)
1603
1604 ! calc the wrapped real coordinates from the wrapped scaled
1605 ! coordinates
1606
1607 d(1) = scaled(1)*Hmat(1,1)
1608 d(2) = scaled(2)*Hmat(2,2)
1609 d(3) = scaled(3)*Hmat(3,3)
1610
1611 endif
1612
1613 endif
1614
1615 r_sq = d(1)*d(1) + d(2)*d(2) + d(3)*d(3)
1616
1617 end subroutine get_interatomic_vector
1618
1619 subroutine zero_work_arrays()
1620
1621 #ifdef IS_MPI
1622
1623 q_Row = 0.0_dp
1624 q_Col = 0.0_dp
1625
1626 q_group_Row = 0.0_dp
1627 q_group_Col = 0.0_dp
1628
1629 eFrame_Row = 0.0_dp
1630 eFrame_Col = 0.0_dp
1631
1632 A_Row = 0.0_dp
1633 A_Col = 0.0_dp
1634
1635 f_Row = 0.0_dp
1636 f_Col = 0.0_dp
1637 f_Temp = 0.0_dp
1638
1639 t_Row = 0.0_dp
1640 t_Col = 0.0_dp
1641 t_Temp = 0.0_dp
1642
1643 pot_Row = 0.0_dp
1644 pot_Col = 0.0_dp
1645 pot_Temp = 0.0_dp
1646
1647 #endif
1648
1649 if (FF_uses_EAM .and. SIM_uses_EAM) then
1650 call clean_EAM()
1651 endif
1652
1653 end subroutine zero_work_arrays
1654
1655 function skipThisPair(atom1, atom2) result(skip_it)
1656 integer, intent(in) :: atom1
1657 integer, intent(in), optional :: atom2
1658 logical :: skip_it
1659 integer :: unique_id_1, unique_id_2
1660 integer :: me_i,me_j
1661 integer :: i
1662
1663 skip_it = .false.
1664
1665 !! there are a number of reasons to skip a pair or a particle
1666 !! mostly we do this to exclude atoms who are involved in short
1667 !! range interactions (bonds, bends, torsions), but we also need
1668 !! to exclude some overcounted interactions that result from
1669 !! the parallel decomposition
1670
1671 #ifdef IS_MPI
1672 !! in MPI, we have to look up the unique IDs for each atom
1673 unique_id_1 = AtomRowToGlobal(atom1)
1674 #else
1675 !! in the normal loop, the atom numbers are unique
1676 unique_id_1 = atom1
1677 #endif
1678
1679 !! We were called with only one atom, so just check the global exclude
1680 !! list for this atom
1681 if (.not. present(atom2)) then
1682 do i = 1, nExcludes_global
1683 if (excludesGlobal(i) == unique_id_1) then
1684 skip_it = .true.
1685 return
1686 end if
1687 end do
1688 return
1689 end if
1690
1691 #ifdef IS_MPI
1692 unique_id_2 = AtomColToGlobal(atom2)
1693 #else
1694 unique_id_2 = atom2
1695 #endif
1696
1697 #ifdef IS_MPI
1698 !! this situation should only arise in MPI simulations
1699 if (unique_id_1 == unique_id_2) then
1700 skip_it = .true.
1701 return
1702 end if
1703
1704 !! this prevents us from doing the pair on multiple processors
1705 if (unique_id_1 < unique_id_2) then
1706 if (mod(unique_id_1 + unique_id_2,2) == 0) then
1707 skip_it = .true.
1708 return
1709 endif
1710 else
1711 if (mod(unique_id_1 + unique_id_2,2) == 1) then
1712 skip_it = .true.
1713 return
1714 endif
1715 endif
1716 #endif
1717
1718 !! the rest of these situations can happen in all simulations:
1719 do i = 1, nExcludes_global
1720 if ((excludesGlobal(i) == unique_id_1) .or. &
1721 (excludesGlobal(i) == unique_id_2)) then
1722 skip_it = .true.
1723 return
1724 endif
1725 enddo
1726
1727 do i = 1, nSkipsForAtom(atom1)
1728 if (skipsForAtom(atom1, i) .eq. unique_id_2) then
1729 skip_it = .true.
1730 return
1731 endif
1732 end do
1733
1734 return
1735 end function skipThisPair
1736
1737 function FF_UsesDirectionalAtoms() result(doesit)
1738 logical :: doesit
1739 doesit = FF_uses_DirectionalAtoms
1740 end function FF_UsesDirectionalAtoms
1741
1742 function FF_RequiresPrepairCalc() result(doesit)
1743 logical :: doesit
1744 doesit = FF_uses_EAM .or. FF_uses_SC &
1745 .or. FF_uses_MEAM
1746 end function FF_RequiresPrepairCalc
1747
1748 #ifdef PROFILE
1749 function getforcetime() result(totalforcetime)
1750 real(kind=dp) :: totalforcetime
1751 totalforcetime = forcetime
1752 end function getforcetime
1753 #endif
1754
1755 !! This cleans componets of force arrays belonging only to fortran
1756
1757 subroutine add_stress_tensor(dpair, fpair, tau)
1758
1759 real( kind = dp ), dimension(3), intent(in) :: dpair, fpair
1760 real( kind = dp ), dimension(9), intent(inout) :: tau
1761
1762 ! because the d vector is the rj - ri vector, and
1763 ! because fx, fy, fz are the force on atom i, we need a
1764 ! negative sign here:
1765
1766 tau(1) = tau(1) - dpair(1) * fpair(1)
1767 tau(2) = tau(2) - dpair(1) * fpair(2)
1768 tau(3) = tau(3) - dpair(1) * fpair(3)
1769 tau(4) = tau(4) - dpair(2) * fpair(1)
1770 tau(5) = tau(5) - dpair(2) * fpair(2)
1771 tau(6) = tau(6) - dpair(2) * fpair(3)
1772 tau(7) = tau(7) - dpair(3) * fpair(1)
1773 tau(8) = tau(8) - dpair(3) * fpair(2)
1774 tau(9) = tau(9) - dpair(3) * fpair(3)
1775
1776 end subroutine add_stress_tensor
1777
1778 end module doForces