ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2284
Committed: Wed Sep 7 19:18:54 2005 UTC (18 years, 10 months ago) by gezelter
File size: 41962 byte(s)
Log Message:
Some bug hunting

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