ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/doForces.F90
Revision: 2310
Committed: Mon Sep 19 23:21:46 2005 UTC (18 years, 9 months ago) by chrisfen
File size: 42105 byte(s)
Log Message:
Fixed bugs in reaction field, now it appears as though it really is working...

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