ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/UseTheForce/doForces.F90
Revision: 2285
Committed: Wed Sep 7 20:46:46 2005 UTC (18 years, 9 months ago) by gezelter
File size: 42015 byte(s)
Log Message:
adding c-side interface to change cutoff Policy

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