ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2533
Committed: Fri Dec 30 23:15:59 2005 UTC (18 years, 6 months ago) by chuckv
File size: 49844 byte(s)
Log Message:
More Sutton-Chen bug fixes.

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