ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 3164
Committed: Thu Jul 12 23:21:00 2007 UTC (17 years, 1 month ago) by chuckv
File size: 58221 byte(s)
Log Message:
More changes to MnM.

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.91 2007-07-12 23:20:00 chuckv Exp $, $Date: 2007-07-12 23:20:00 $, $Name: not supported by cvs2svn $, $Revision: 1.91 $
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, &
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 logical ( kind = 2) :: do_pot_c, do_stress_c
816 logical :: do_pot
817 logical :: do_stress
818 logical :: in_switching_region
819 #ifdef IS_MPI
820 real( kind = DP ), dimension(LR_POT_TYPES) :: pot_local
821 integer :: nAtomsInRow
822 integer :: nAtomsInCol
823 integer :: nprocs
824 integer :: nGroupsInRow
825 integer :: nGroupsInCol
826 #endif
827 integer :: natoms
828 logical :: update_nlist
829 integer :: i, j, jstart, jend, jnab
830 integer :: istart, iend
831 integer :: ia, jb, atom1, atom2
832 integer :: nlist
833 real( kind = DP ) :: ratmsq, rgrpsq, rgrp, rag, vpair, vij
834 real( kind = DP ) :: sw, dswdr, swderiv, mf
835 real( kind = DP ) :: rVal
836 real(kind=dp),dimension(3) :: d_atm, d_grp, fpair, fij, fg, dag
837 real(kind=dp) :: rfpot, mu_i
838 real(kind=dp):: rCut
839 integer :: me_i, me_j, n_in_i, n_in_j
840 logical :: is_dp_i
841 integer :: neighborListSize
842 integer :: listerror, error
843 integer :: localError
844 integer :: propPack_i, propPack_j
845 integer :: loopStart, loopEnd, loop
846 integer :: iHash
847 integer :: i1
848
849 !! the variables for the box dipole moment
850 #ifdef IS_MPI
851 integer :: pChgCount_local
852 integer :: nChgCount_local
853 real(kind=dp) :: pChg_local
854 real(kind=dp) :: nChg_local
855 real(kind=dp), dimension(3) :: pChgPos_local
856 real(kind=dp), dimension(3) :: nChgPos_local
857 real(kind=dp), dimension(3) :: dipVec_local
858 #endif
859 integer :: pChgCount
860 integer :: nChgCount
861 real(kind=dp) :: pChg
862 real(kind=dp) :: nChg
863 real(kind=dp) :: chg_value
864 real(kind=dp), dimension(3) :: pChgPos
865 real(kind=dp), dimension(3) :: nChgPos
866 real(kind=dp), dimension(3) :: dipVec
867 real(kind=dp), dimension(3) :: chgVec
868
869 !! initialize box dipole variables
870 if (do_box_dipole) then
871 #ifdef IS_MPI
872 pChg_local = 0.0_dp
873 nChg_local = 0.0_dp
874 pChgCount_local = 0
875 nChgCount_local = 0
876 do i=1, 3
877 pChgPos_local = 0.0_dp
878 nChgPos_local = 0.0_dp
879 dipVec_local = 0.0_dp
880 enddo
881 #endif
882 pChg = 0.0_dp
883 nChg = 0.0_dp
884 pChgCount = 0
885 nChgCount = 0
886 chg_value = 0.0_dp
887
888 do i=1, 3
889 pChgPos(i) = 0.0_dp
890 nChgPos(i) = 0.0_dp
891 dipVec(i) = 0.0_dp
892 chgVec(i) = 0.0_dp
893 boxDipole(i) = 0.0_dp
894 enddo
895 endif
896
897 !! initialize local variables
898
899 #ifdef IS_MPI
900 pot_local = 0.0_dp
901 nAtomsInRow = getNatomsInRow(plan_atom_row)
902 nAtomsInCol = getNatomsInCol(plan_atom_col)
903 nGroupsInRow = getNgroupsInRow(plan_group_row)
904 nGroupsInCol = getNgroupsInCol(plan_group_col)
905 #else
906 natoms = nlocal
907 #endif
908
909 call doReadyCheck(localError)
910 if ( localError .ne. 0 ) then
911 call handleError("do_force_loop", "Not Initialized")
912 error = -1
913 return
914 end if
915 call zero_work_arrays()
916
917 do_pot = do_pot_c
918 do_stress = do_stress_c
919
920 ! Gather all information needed by all force loops:
921
922 #ifdef IS_MPI
923
924 call gather(q, q_Row, plan_atom_row_3d)
925 call gather(q, q_Col, plan_atom_col_3d)
926
927 call gather(q_group, q_group_Row, plan_group_row_3d)
928 call gather(q_group, q_group_Col, plan_group_col_3d)
929
930 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
931 call gather(eFrame, eFrame_Row, plan_atom_row_rotation)
932 call gather(eFrame, eFrame_Col, plan_atom_col_rotation)
933
934 call gather(A, A_Row, plan_atom_row_rotation)
935 call gather(A, A_Col, plan_atom_col_rotation)
936 endif
937
938 #endif
939
940 !! Begin force loop timing:
941 #ifdef PROFILE
942 call cpu_time(forceTimeInitial)
943 nloops = nloops + 1
944 #endif
945
946 loopEnd = PAIR_LOOP
947 if (FF_RequiresPrepairCalc() .and. SIM_requires_prepair_calc) then
948 loopStart = PREPAIR_LOOP
949 else
950 loopStart = PAIR_LOOP
951 endif
952
953 do loop = loopStart, loopEnd
954
955 ! See if we need to update neighbor lists
956 ! (but only on the first time through):
957 if (loop .eq. loopStart) then
958 #ifdef IS_MPI
959 call checkNeighborList(nGroupsInRow, q_group_row, skinThickness, &
960 update_nlist)
961 #else
962 call checkNeighborList(nGroups, q_group, skinThickness, &
963 update_nlist)
964 #endif
965 endif
966
967 if (update_nlist) then
968 !! save current configuration and construct neighbor list
969 #ifdef IS_MPI
970 call saveNeighborList(nGroupsInRow, q_group_row)
971 #else
972 call saveNeighborList(nGroups, q_group)
973 #endif
974 neighborListSize = size(list)
975 nlist = 0
976 endif
977
978 istart = 1
979 #ifdef IS_MPI
980 iend = nGroupsInRow
981 #else
982 iend = nGroups - 1
983 #endif
984 outer: do i = istart, iend
985
986 if (update_nlist) point(i) = nlist + 1
987
988 n_in_i = groupStartRow(i+1) - groupStartRow(i)
989
990 if (update_nlist) then
991 #ifdef IS_MPI
992 jstart = 1
993 jend = nGroupsInCol
994 #else
995 jstart = i+1
996 jend = nGroups
997 #endif
998 else
999 jstart = point(i)
1000 jend = point(i+1) - 1
1001 ! make sure group i has neighbors
1002 if (jstart .gt. jend) cycle outer
1003 endif
1004
1005 do jnab = jstart, jend
1006 if (update_nlist) then
1007 j = jnab
1008 else
1009 j = list(jnab)
1010 endif
1011
1012 #ifdef IS_MPI
1013 me_j = atid_col(j)
1014 call get_interatomic_vector(q_group_Row(:,i), &
1015 q_group_Col(:,j), d_grp, rgrpsq)
1016 #else
1017 me_j = atid(j)
1018 call get_interatomic_vector(q_group(:,i), &
1019 q_group(:,j), d_grp, rgrpsq)
1020 #endif
1021
1022 if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rListsq) then
1023 if (update_nlist) then
1024 nlist = nlist + 1
1025
1026 if (nlist > neighborListSize) then
1027 #ifdef IS_MPI
1028 call expandNeighborList(nGroupsInRow, listerror)
1029 #else
1030 call expandNeighborList(nGroups, listerror)
1031 #endif
1032 if (listerror /= 0) then
1033 error = -1
1034 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
1035 return
1036 end if
1037 neighborListSize = size(list)
1038 endif
1039
1040 list(nlist) = j
1041 endif
1042
1043 if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCutsq) then
1044
1045 rCut = gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCut
1046 if (loop .eq. PAIR_LOOP) then
1047 vij = 0.0_dp
1048 fij(1) = 0.0_dp
1049 fij(2) = 0.0_dp
1050 fij(3) = 0.0_dp
1051 endif
1052
1053 call get_switch(rgrpsq, sw, dswdr,rgrp, in_switching_region)
1054
1055 n_in_j = groupStartCol(j+1) - groupStartCol(j)
1056
1057 do ia = groupStartRow(i), groupStartRow(i+1)-1
1058
1059 atom1 = groupListRow(ia)
1060
1061 inner: do jb = groupStartCol(j), groupStartCol(j+1)-1
1062
1063 atom2 = groupListCol(jb)
1064
1065 if (skipThisPair(atom1, atom2)) cycle inner
1066
1067 if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
1068 d_atm(1) = d_grp(1)
1069 d_atm(2) = d_grp(2)
1070 d_atm(3) = d_grp(3)
1071 ratmsq = rgrpsq
1072 else
1073 #ifdef IS_MPI
1074 call get_interatomic_vector(q_Row(:,atom1), &
1075 q_Col(:,atom2), d_atm, ratmsq)
1076 #else
1077 call get_interatomic_vector(q(:,atom1), &
1078 q(:,atom2), d_atm, ratmsq)
1079 #endif
1080 endif
1081
1082 if (loop .eq. PREPAIR_LOOP) then
1083 #ifdef IS_MPI
1084 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
1085 rgrpsq, d_grp, rCut, do_pot, do_stress, &
1086 eFrame, A, f, t, pot_local)
1087 #else
1088 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
1089 rgrpsq, d_grp, rCut, do_pot, do_stress, &
1090 eFrame, A, f, t, pot)
1091 #endif
1092 else
1093 #ifdef IS_MPI
1094 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
1095 do_pot, eFrame, A, f, t, pot_local, vpair, &
1096 fpair, d_grp, rgrp, rCut)
1097 #else
1098 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
1099 do_pot, eFrame, A, f, t, pot, vpair, fpair, &
1100 d_grp, rgrp, rCut)
1101 #endif
1102 vij = vij + vpair
1103 fij(1) = fij(1) + fpair(1)
1104 fij(2) = fij(2) + fpair(2)
1105 fij(3) = fij(3) + fpair(3)
1106 if (do_stress) then
1107 call add_stress_tensor(d_atm, fpair, tau)
1108 endif
1109 endif
1110 enddo inner
1111 enddo
1112
1113 if (loop .eq. PAIR_LOOP) then
1114 if (in_switching_region) then
1115 swderiv = vij*dswdr/rgrp
1116 fg = swderiv*d_grp
1117
1118 fij(1) = fij(1) + fg(1)
1119 fij(2) = fij(2) + fg(2)
1120 fij(3) = fij(3) + fg(3)
1121
1122 if (do_stress .and. (n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
1123 call add_stress_tensor(d_atm, fg, tau)
1124 endif
1125
1126 do ia=groupStartRow(i), groupStartRow(i+1)-1
1127 atom1=groupListRow(ia)
1128 mf = mfactRow(atom1)
1129 ! fg is the force on atom ia due to cutoff group's
1130 ! presence in switching region
1131 fg = swderiv*d_grp*mf
1132 #ifdef IS_MPI
1133 f_Row(1,atom1) = f_Row(1,atom1) + fg(1)
1134 f_Row(2,atom1) = f_Row(2,atom1) + fg(2)
1135 f_Row(3,atom1) = f_Row(3,atom1) + fg(3)
1136 #else
1137 f(1,atom1) = f(1,atom1) + fg(1)
1138 f(2,atom1) = f(2,atom1) + fg(2)
1139 f(3,atom1) = f(3,atom1) + fg(3)
1140 #endif
1141 if (n_in_i .gt. 1) then
1142 if (do_stress.and.SIM_uses_AtomicVirial) then
1143 ! find the distance between the atom and the center of
1144 ! the cutoff group:
1145 #ifdef IS_MPI
1146 call get_interatomic_vector(q_Row(:,atom1), &
1147 q_group_Row(:,i), dag, rag)
1148 #else
1149 call get_interatomic_vector(q(:,atom1), &
1150 q_group(:,i), dag, rag)
1151 #endif
1152 call add_stress_tensor(dag,fg,tau)
1153 endif
1154 endif
1155 enddo
1156
1157 do jb=groupStartCol(j), groupStartCol(j+1)-1
1158 atom2=groupListCol(jb)
1159 mf = mfactCol(atom2)
1160 ! fg is the force on atom jb due to cutoff group's
1161 ! presence in switching region
1162 fg = -swderiv*d_grp*mf
1163 #ifdef IS_MPI
1164 f_Col(1,atom2) = f_Col(1,atom2) + fg(1)
1165 f_Col(2,atom2) = f_Col(2,atom2) + fg(2)
1166 f_Col(3,atom2) = f_Col(3,atom2) + fg(3)
1167 #else
1168 f(1,atom2) = f(1,atom2) + fg(1)
1169 f(2,atom2) = f(2,atom2) + fg(2)
1170 f(3,atom2) = f(3,atom2) + fg(3)
1171 #endif
1172 if (n_in_j .gt. 1) then
1173 if (do_stress.and.SIM_uses_AtomicVirial) then
1174 ! find the distance between the atom and the center of
1175 ! the cutoff group:
1176 #ifdef IS_MPI
1177 call get_interatomic_vector(q_Col(:,atom2), &
1178 q_group_Col(:,j), dag, rag)
1179 #else
1180 call get_interatomic_vector(q(:,atom2), &
1181 q_group(:,j), dag, rag)
1182 #endif
1183 call add_stress_tensor(dag,fg,tau)
1184 endif
1185 endif
1186 enddo
1187 endif
1188 endif
1189 endif
1190 endif
1191 enddo
1192
1193 enddo outer
1194
1195 if (update_nlist) then
1196 #ifdef IS_MPI
1197 point(nGroupsInRow + 1) = nlist + 1
1198 #else
1199 point(nGroups) = nlist + 1
1200 #endif
1201 if (loop .eq. PREPAIR_LOOP) then
1202 ! we just did the neighbor list update on the first
1203 ! pass, so we don't need to do it
1204 ! again on the second pass
1205 update_nlist = .false.
1206 endif
1207 endif
1208
1209 if (loop .eq. PREPAIR_LOOP) then
1210 #ifdef IS_MPI
1211 call do_preforce(nlocal, pot_local)
1212 #else
1213 call do_preforce(nlocal, pot)
1214 #endif
1215 endif
1216
1217 enddo
1218
1219 !! Do timing
1220 #ifdef PROFILE
1221 call cpu_time(forceTimeFinal)
1222 forceTime = forceTime + forceTimeFinal - forceTimeInitial
1223 #endif
1224
1225 #ifdef IS_MPI
1226 !!distribute forces
1227
1228 f_temp = 0.0_dp
1229 call scatter(f_Row,f_temp,plan_atom_row_3d)
1230 do i = 1,nlocal
1231 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1232 end do
1233
1234 f_temp = 0.0_dp
1235 call scatter(f_Col,f_temp,plan_atom_col_3d)
1236 do i = 1,nlocal
1237 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1238 end do
1239
1240 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
1241 t_temp = 0.0_dp
1242 call scatter(t_Row,t_temp,plan_atom_row_3d)
1243 do i = 1,nlocal
1244 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1245 end do
1246 t_temp = 0.0_dp
1247 call scatter(t_Col,t_temp,plan_atom_col_3d)
1248
1249 do i = 1,nlocal
1250 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1251 end do
1252 endif
1253
1254 if (do_pot) then
1255 ! scatter/gather pot_row into the members of my column
1256 do i = 1,LR_POT_TYPES
1257 call scatter(pot_Row(i,:), pot_Temp(i,:), plan_atom_row)
1258 end do
1259 ! scatter/gather pot_local into all other procs
1260 ! add resultant to get total pot
1261 do i = 1, nlocal
1262 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES) &
1263 + pot_Temp(1:LR_POT_TYPES,i)
1264 enddo
1265
1266 pot_Temp = 0.0_DP
1267 do i = 1,LR_POT_TYPES
1268 call scatter(pot_Col(i,:), pot_Temp(i,:), plan_atom_col)
1269 end do
1270 do i = 1, nlocal
1271 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES)&
1272 + pot_Temp(1:LR_POT_TYPES,i)
1273 enddo
1274
1275 endif
1276 #endif
1277
1278 if (SIM_requires_postpair_calc) then
1279 do i = 1, nlocal
1280
1281 ! we loop only over the local atoms, so we don't need row and column
1282 ! lookups for the types
1283
1284 me_i = atid(i)
1285
1286 ! is the atom electrostatic? See if it would have an
1287 ! electrostatic interaction with itself
1288 iHash = InteractionHash(me_i,me_i)
1289
1290 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1291 #ifdef IS_MPI
1292 call self_self(i, eFrame, pot_local(ELECTROSTATIC_POT), &
1293 t, do_pot)
1294 #else
1295 call self_self(i, eFrame, pot(ELECTROSTATIC_POT), &
1296 t, do_pot)
1297 #endif
1298 endif
1299
1300
1301 if (electrostaticSummationMethod.eq.REACTION_FIELD) then
1302
1303 ! loop over the excludes to accumulate RF stuff we've
1304 ! left out of the normal pair loop
1305
1306 do i1 = 1, nSkipsForAtom(i)
1307 j = skipsForAtom(i, i1)
1308
1309 ! prevent overcounting of the skips
1310 if (i.lt.j) then
1311 call get_interatomic_vector(q(:,i), q(:,j), d_atm, ratmsq)
1312 rVal = sqrt(ratmsq)
1313 call get_switch(ratmsq, sw, dswdr, rVal,in_switching_region)
1314 #ifdef IS_MPI
1315 call rf_self_excludes(i, j, sw, eFrame, d_atm, rVal, &
1316 vpair, pot_local(ELECTROSTATIC_POT), f, t, do_pot)
1317 #else
1318 call rf_self_excludes(i, j, sw, eFrame, d_atm, rVal, &
1319 vpair, pot(ELECTROSTATIC_POT), f, t, do_pot)
1320 #endif
1321 endif
1322 enddo
1323 endif
1324
1325 if (do_box_dipole) then
1326 #ifdef IS_MPI
1327 call accumulate_box_dipole(i, eFrame, q(:,i), pChg_local, &
1328 nChg_local, pChgPos_local, nChgPos_local, dipVec_local, &
1329 pChgCount_local, nChgCount_local)
1330 #else
1331 call accumulate_box_dipole(i, eFrame, q(:,i), pChg, nChg, &
1332 pChgPos, nChgPos, dipVec, pChgCount, nChgCount)
1333 #endif
1334 endif
1335 enddo
1336 endif
1337
1338 #ifdef IS_MPI
1339 if (do_pot) then
1340 #ifdef SINGLE_PRECISION
1341 call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_real,mpi_sum, &
1342 mpi_comm_world,mpi_err)
1343 #else
1344 call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_double_precision, &
1345 mpi_sum, mpi_comm_world,mpi_err)
1346 #endif
1347 endif
1348
1349 if (do_box_dipole) then
1350
1351 #ifdef SINGLE_PRECISION
1352 call mpi_allreduce(pChg_local, pChg, 1, mpi_real, mpi_sum, &
1353 mpi_comm_world, mpi_err)
1354 call mpi_allreduce(nChg_local, nChg, 1, mpi_real, mpi_sum, &
1355 mpi_comm_world, mpi_err)
1356 call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer, mpi_sum,&
1357 mpi_comm_world, mpi_err)
1358 call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer, mpi_sum,&
1359 mpi_comm_world, mpi_err)
1360 call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_real, mpi_sum, &
1361 mpi_comm_world, mpi_err)
1362 call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_real, mpi_sum, &
1363 mpi_comm_world, mpi_err)
1364 call mpi_allreduce(dipVec_local, dipVec, 3, mpi_real, mpi_sum, &
1365 mpi_comm_world, mpi_err)
1366 #else
1367 call mpi_allreduce(pChg_local, pChg, 1, mpi_double_precision, mpi_sum, &
1368 mpi_comm_world, mpi_err)
1369 call mpi_allreduce(nChg_local, nChg, 1, mpi_double_precision, mpi_sum, &
1370 mpi_comm_world, mpi_err)
1371 call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer,&
1372 mpi_sum, mpi_comm_world, mpi_err)
1373 call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer,&
1374 mpi_sum, mpi_comm_world, mpi_err)
1375 call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_double_precision, &
1376 mpi_sum, mpi_comm_world, mpi_err)
1377 call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_double_precision, &
1378 mpi_sum, mpi_comm_world, mpi_err)
1379 call mpi_allreduce(dipVec_local, dipVec, 3, mpi_double_precision, &
1380 mpi_sum, mpi_comm_world, mpi_err)
1381 #endif
1382
1383 endif
1384
1385 #endif
1386
1387 if (do_box_dipole) then
1388 ! first load the accumulated dipole moment (if dipoles were present)
1389 boxDipole(1) = dipVec(1)
1390 boxDipole(2) = dipVec(2)
1391 boxDipole(3) = dipVec(3)
1392
1393 ! now include the dipole moment due to charges
1394 ! use the lesser of the positive and negative charge totals
1395 if (nChg .le. pChg) then
1396 chg_value = nChg
1397 else
1398 chg_value = pChg
1399 endif
1400
1401 ! find the average positions
1402 if (pChgCount .gt. 0 .and. nChgCount .gt. 0) then
1403 pChgPos = pChgPos / pChgCount
1404 nChgPos = nChgPos / nChgCount
1405 endif
1406
1407 ! dipole is from the negative to the positive (physics notation)
1408 chgVec(1) = pChgPos(1) - nChgPos(1)
1409 chgVec(2) = pChgPos(2) - nChgPos(2)
1410 chgVec(3) = pChgPos(3) - nChgPos(3)
1411
1412 boxDipole(1) = boxDipole(1) + chgVec(1) * chg_value
1413 boxDipole(2) = boxDipole(2) + chgVec(2) * chg_value
1414 boxDipole(3) = boxDipole(3) + chgVec(3) * chg_value
1415
1416 endif
1417
1418 end subroutine do_force_loop
1419
1420 subroutine do_pair(i, j, rijsq, d, sw, do_pot, &
1421 eFrame, A, f, t, pot, vpair, fpair, d_grp, r_grp, rCut)
1422
1423 real( kind = dp ) :: vpair, sw
1424 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1425 real( kind = dp ), dimension(3) :: fpair
1426 real( kind = dp ), dimension(nLocal) :: mfact
1427 real( kind = dp ), dimension(9,nLocal) :: eFrame
1428 real( kind = dp ), dimension(9,nLocal) :: A
1429 real( kind = dp ), dimension(3,nLocal) :: f
1430 real( kind = dp ), dimension(3,nLocal) :: t
1431
1432 logical, intent(inout) :: do_pot
1433 integer, intent(in) :: i, j
1434 real ( kind = dp ), intent(inout) :: rijsq
1435 real ( kind = dp ), intent(inout) :: r_grp
1436 real ( kind = dp ), intent(inout) :: d(3)
1437 real ( kind = dp ), intent(inout) :: d_grp(3)
1438 real ( kind = dp ), intent(inout) :: rCut
1439 real ( kind = dp ) :: r
1440 real ( kind = dp ) :: a_k, b_k, c_k, d_k, dx
1441 integer :: me_i, me_j
1442 integer :: k
1443
1444 integer :: iHash
1445
1446 r = sqrt(rijsq)
1447
1448 vpair = 0.0_dp
1449 fpair(1:3) = 0.0_dp
1450
1451 #ifdef IS_MPI
1452 me_i = atid_row(i)
1453 me_j = atid_col(j)
1454 #else
1455 me_i = atid(i)
1456 me_j = atid(j)
1457 #endif
1458
1459 iHash = InteractionHash(me_i, me_j)
1460
1461 if ( iand(iHash, LJ_PAIR).ne.0 ) then
1462 call do_lj_pair(i, j, d, r, rijsq, rcut, sw, vpair, fpair, &
1463 pot(VDW_POT), f, do_pot)
1464 endif
1465
1466 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1467 call doElectrostaticPair(i, j, d, r, rijsq, rcut, sw, vpair, fpair, &
1468 pot(ELECTROSTATIC_POT), eFrame, f, t, do_pot)
1469 endif
1470
1471 if ( iand(iHash, STICKY_PAIR).ne.0 ) then
1472 call do_sticky_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1473 pot(HB_POT), A, f, t, do_pot)
1474 endif
1475
1476 if ( iand(iHash, STICKYPOWER_PAIR).ne.0 ) then
1477 call do_sticky_power_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1478 pot(HB_POT), A, f, t, do_pot)
1479 endif
1480
1481 if ( iand(iHash, GAYBERNE_PAIR).ne.0 ) then
1482 call do_gb_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1483 pot(VDW_POT), A, f, t, do_pot)
1484 endif
1485
1486 if ( iand(iHash, GAYBERNE_LJ).ne.0 ) then
1487 call do_gb_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1488 pot(VDW_POT), A, f, t, do_pot)
1489 endif
1490
1491 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1492 call do_eam_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1493 pot(METALLIC_POT), f, do_pot)
1494 endif
1495
1496 if ( iand(iHash, SHAPE_PAIR).ne.0 ) then
1497 call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1498 pot(VDW_POT), A, f, t, do_pot)
1499 endif
1500
1501 if ( iand(iHash, SHAPE_LJ).ne.0 ) then
1502 call do_shape_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, SC_PAIR).ne.0 ) then
1507 call do_SC_pair(i, j, d, r, rijsq, rcut, sw, vpair, fpair, &
1508 pot(METALLIC_POT), f, do_pot)
1509 endif
1510
1511 end subroutine do_pair
1512
1513 subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, rCut, &
1514 do_pot, do_stress, eFrame, A, f, t, pot)
1515
1516 real( kind = dp ) :: sw
1517 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1518 real( kind = dp ), dimension(9,nLocal) :: eFrame
1519 real (kind=dp), dimension(9,nLocal) :: A
1520 real (kind=dp), dimension(3,nLocal) :: f
1521 real (kind=dp), dimension(3,nLocal) :: t
1522
1523 logical, intent(inout) :: do_pot, do_stress
1524 integer, intent(in) :: i, j
1525 real ( kind = dp ), intent(inout) :: rijsq, rcijsq, rCut
1526 real ( kind = dp ) :: r, rc
1527 real ( kind = dp ), intent(inout) :: d(3), dc(3)
1528
1529 integer :: me_i, me_j, iHash
1530
1531 r = sqrt(rijsq)
1532
1533 #ifdef IS_MPI
1534 me_i = atid_row(i)
1535 me_j = atid_col(j)
1536 #else
1537 me_i = atid(i)
1538 me_j = atid(j)
1539 #endif
1540
1541 iHash = InteractionHash(me_i, me_j)
1542
1543 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1544 call calc_EAM_prepair_rho(i, j, d, r, rijsq)
1545 endif
1546
1547 if ( iand(iHash, SC_PAIR).ne.0 ) then
1548 call calc_SC_prepair_rho(i, j, d, r, rijsq, rcut )
1549 endif
1550
1551 end subroutine do_prepair
1552
1553
1554 subroutine do_preforce(nlocal,pot)
1555 integer :: nlocal
1556 real( kind = dp ),dimension(LR_POT_TYPES) :: pot
1557
1558 if (FF_uses_EAM .and. SIM_uses_EAM) then
1559 call calc_EAM_preforce_Frho(nlocal,pot(METALLIC_POT))
1560 endif
1561 if (FF_uses_SC .and. SIM_uses_SC) then
1562 call calc_SC_preforce_Frho(nlocal,pot(METALLIC_POT))
1563 endif
1564 end subroutine do_preforce
1565
1566
1567 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
1568
1569 real (kind = dp), dimension(3) :: q_i
1570 real (kind = dp), dimension(3) :: q_j
1571 real ( kind = dp ), intent(out) :: r_sq
1572 real( kind = dp ) :: d(3), scaled(3)
1573 integer i
1574
1575 d(1) = q_j(1) - q_i(1)
1576 d(2) = q_j(2) - q_i(2)
1577 d(3) = q_j(3) - q_i(3)
1578
1579 ! Wrap back into periodic box if necessary
1580 if ( SIM_uses_PBC ) then
1581
1582 if( .not.boxIsOrthorhombic ) then
1583 ! calc the scaled coordinates.
1584 ! scaled = matmul(HmatInv, d)
1585
1586 scaled(1) = HmatInv(1,1)*d(1) + HmatInv(1,2)*d(2) + HmatInv(1,3)*d(3)
1587 scaled(2) = HmatInv(2,1)*d(1) + HmatInv(2,2)*d(2) + HmatInv(2,3)*d(3)
1588 scaled(3) = HmatInv(3,1)*d(1) + HmatInv(3,2)*d(2) + HmatInv(3,3)*d(3)
1589
1590 ! wrap the scaled coordinates
1591
1592 scaled(1) = scaled(1) - anint(scaled(1), kind=dp)
1593 scaled(2) = scaled(2) - anint(scaled(2), kind=dp)
1594 scaled(3) = scaled(3) - anint(scaled(3), kind=dp)
1595
1596 ! calc the wrapped real coordinates from the wrapped scaled
1597 ! coordinates
1598 ! d = matmul(Hmat,scaled)
1599 d(1)= Hmat(1,1)*scaled(1) + Hmat(1,2)*scaled(2) + Hmat(1,3)*scaled(3)
1600 d(2)= Hmat(2,1)*scaled(1) + Hmat(2,2)*scaled(2) + Hmat(2,3)*scaled(3)
1601 d(3)= Hmat(3,1)*scaled(1) + Hmat(3,2)*scaled(2) + Hmat(3,3)*scaled(3)
1602
1603 else
1604 ! calc the scaled coordinates.
1605
1606 scaled(1) = d(1) * HmatInv(1,1)
1607 scaled(2) = d(2) * HmatInv(2,2)
1608 scaled(3) = d(3) * HmatInv(3,3)
1609
1610 ! wrap the scaled coordinates
1611
1612 scaled(1) = scaled(1) - anint(scaled(1), kind=dp)
1613 scaled(2) = scaled(2) - anint(scaled(2), kind=dp)
1614 scaled(3) = scaled(3) - anint(scaled(3), kind=dp)
1615
1616 ! calc the wrapped real coordinates from the wrapped scaled
1617 ! coordinates
1618
1619 d(1) = scaled(1)*Hmat(1,1)
1620 d(2) = scaled(2)*Hmat(2,2)
1621 d(3) = scaled(3)*Hmat(3,3)
1622
1623 endif
1624
1625 endif
1626
1627 r_sq = d(1)*d(1) + d(2)*d(2) + d(3)*d(3)
1628
1629 end subroutine get_interatomic_vector
1630
1631 subroutine zero_work_arrays()
1632
1633 #ifdef IS_MPI
1634
1635 q_Row = 0.0_dp
1636 q_Col = 0.0_dp
1637
1638 q_group_Row = 0.0_dp
1639 q_group_Col = 0.0_dp
1640
1641 eFrame_Row = 0.0_dp
1642 eFrame_Col = 0.0_dp
1643
1644 A_Row = 0.0_dp
1645 A_Col = 0.0_dp
1646
1647 f_Row = 0.0_dp
1648 f_Col = 0.0_dp
1649 f_Temp = 0.0_dp
1650
1651 t_Row = 0.0_dp
1652 t_Col = 0.0_dp
1653 t_Temp = 0.0_dp
1654
1655 pot_Row = 0.0_dp
1656 pot_Col = 0.0_dp
1657 pot_Temp = 0.0_dp
1658
1659 #endif
1660
1661 if (FF_uses_EAM .and. SIM_uses_EAM) then
1662 call clean_EAM()
1663 endif
1664
1665 end subroutine zero_work_arrays
1666
1667 function skipThisPair(atom1, atom2) result(skip_it)
1668 integer, intent(in) :: atom1
1669 integer, intent(in), optional :: atom2
1670 logical :: skip_it
1671 integer :: unique_id_1, unique_id_2
1672 integer :: me_i,me_j
1673 integer :: i
1674
1675 skip_it = .false.
1676
1677 !! there are a number of reasons to skip a pair or a particle
1678 !! mostly we do this to exclude atoms who are involved in short
1679 !! range interactions (bonds, bends, torsions), but we also need
1680 !! to exclude some overcounted interactions that result from
1681 !! the parallel decomposition
1682
1683 #ifdef IS_MPI
1684 !! in MPI, we have to look up the unique IDs for each atom
1685 unique_id_1 = AtomRowToGlobal(atom1)
1686 #else
1687 !! in the normal loop, the atom numbers are unique
1688 unique_id_1 = atom1
1689 #endif
1690
1691 !! We were called with only one atom, so just check the global exclude
1692 !! list for this atom
1693 if (.not. present(atom2)) then
1694 do i = 1, nExcludes_global
1695 if (excludesGlobal(i) == unique_id_1) then
1696 skip_it = .true.
1697 return
1698 end if
1699 end do
1700 return
1701 end if
1702
1703 #ifdef IS_MPI
1704 unique_id_2 = AtomColToGlobal(atom2)
1705 #else
1706 unique_id_2 = atom2
1707 #endif
1708
1709 #ifdef IS_MPI
1710 !! this situation should only arise in MPI simulations
1711 if (unique_id_1 == unique_id_2) then
1712 skip_it = .true.
1713 return
1714 end if
1715
1716 !! this prevents us from doing the pair on multiple processors
1717 if (unique_id_1 < unique_id_2) then
1718 if (mod(unique_id_1 + unique_id_2,2) == 0) then
1719 skip_it = .true.
1720 return
1721 endif
1722 else
1723 if (mod(unique_id_1 + unique_id_2,2) == 1) then
1724 skip_it = .true.
1725 return
1726 endif
1727 endif
1728 #endif
1729
1730 !! the rest of these situations can happen in all simulations:
1731 do i = 1, nExcludes_global
1732 if ((excludesGlobal(i) == unique_id_1) .or. &
1733 (excludesGlobal(i) == unique_id_2)) then
1734 skip_it = .true.
1735 return
1736 endif
1737 enddo
1738
1739 do i = 1, nSkipsForAtom(atom1)
1740 if (skipsForAtom(atom1, i) .eq. unique_id_2) then
1741 skip_it = .true.
1742 return
1743 endif
1744 end do
1745
1746 return
1747 end function skipThisPair
1748
1749 function FF_UsesDirectionalAtoms() result(doesit)
1750 logical :: doesit
1751 doesit = FF_uses_DirectionalAtoms
1752 end function FF_UsesDirectionalAtoms
1753
1754 function FF_RequiresPrepairCalc() result(doesit)
1755 logical :: doesit
1756 doesit = FF_uses_EAM .or. FF_uses_SC
1757 end function FF_RequiresPrepairCalc
1758
1759 #ifdef PROFILE
1760 function getforcetime() result(totalforcetime)
1761 real(kind=dp) :: totalforcetime
1762 totalforcetime = forcetime
1763 end function getforcetime
1764 #endif
1765
1766 !! This cleans componets of force arrays belonging only to fortran
1767
1768 subroutine add_stress_tensor(dpair, fpair, tau)
1769
1770 real( kind = dp ), dimension(3), intent(in) :: dpair, fpair
1771 real( kind = dp ), dimension(9), intent(inout) :: tau
1772
1773 ! because the d vector is the rj - ri vector, and
1774 ! because fx, fy, fz are the force on atom i, we need a
1775 ! negative sign here:
1776
1777 tau(1) = tau(1) - dpair(1) * fpair(1)
1778 tau(2) = tau(2) - dpair(1) * fpair(2)
1779 tau(3) = tau(3) - dpair(1) * fpair(3)
1780 tau(4) = tau(4) - dpair(2) * fpair(1)
1781 tau(5) = tau(5) - dpair(2) * fpair(2)
1782 tau(6) = tau(6) - dpair(2) * fpair(3)
1783 tau(7) = tau(7) - dpair(3) * fpair(1)
1784 tau(8) = tau(8) - dpair(3) * fpair(2)
1785 tau(9) = tau(9) - dpair(3) * fpair(3)
1786
1787 end subroutine add_stress_tensor
1788
1789 end module doForces