ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/doForces.F90
Revision: 2290
Committed: Wed Sep 14 20:28:05 2005 UTC (18 years, 9 months ago) by gezelter
File size: 42224 byte(s)
Log Message:
fix to put back calculation of r in do_prepair

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