ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2323
Committed: Fri Sep 23 20:31:02 2005 UTC (18 years, 9 months ago) by chuckv
File size: 42760 byte(s)
Log Message:
Fixed "dum-dum" where we ignore the skin thickness and hardcode listSkin to be 1.0. We now get listskin from skin. This will get fixed to where we can manually set skin thickness.

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