ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2917
Committed: Mon Jul 3 13:18:43 2006 UTC (18 years ago) by chrisfen
File size: 56113 byte(s)
Log Message:
Added simulation box dipole moment accumulation for the purposes of calculating dielectric constants

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