ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2282
Committed: Tue Sep 6 17:32:42 2005 UTC (18 years, 10 months ago) by chuckv
File size: 41959 byte(s)
Log Message:
Added allocation for gtypeCutoffmap etc..

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