ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2281
Committed: Thu Sep 1 22:56:20 2005 UTC (18 years, 10 months ago) by gezelter
File size: 41743 byte(s)
Log Message:
initialized atomTypeMaxCutoff(i) to zero

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.35 2005-09-01 22:56:20 gezelter Exp $, $Date: 2005-09-01 22:56:20 $, $Name: not supported by cvs2svn $, $Revision: 1.35 $
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
334 !! first we do a single loop over the cutoff groups to find the
335 !! largest cutoff for any atypes present in this group. We also
336 !! create gtypes at this point.
337
338 tol = 1.0d-6
339
340 do i = istart, iend
341 n_in_i = groupStartRow(i+1) - groupStartRow(i)
342 groupMaxCutoff(i) = 0.0_dp
343 do ia = groupStartRow(i), groupStartRow(i+1)-1
344 atom1 = groupListRow(ia)
345 #ifdef IS_MPI
346 me_i = atid_row(atom1)
347 #else
348 me_i = atid(atom1)
349 #endif
350 if (atypeMaxCutoff(me_i).gt.groupMaxCutoff(i)) then
351 groupMaxCutoff(i)=atypeMaxCutoff(me_i)
352 endif
353 enddo
354 if (nGroupTypes.eq.0) then
355 nGroupTypes = nGroupTypes + 1
356 gtypeMaxCutoff(nGroupTypes) = groupMaxCutoff(i)
357 groupToGtype(i) = nGroupTypes
358 else
359 do g = 1, nGroupTypes
360 if ( abs(groupMaxCutoff(i) - gtypeMaxCutoff(g)).gt.tol) then
361 nGroupTypes = nGroupTypes + 1
362 gtypeMaxCutoff(nGroupTypes) = groupMaxCutoff(i)
363 groupToGtype(i) = nGroupTypes
364 else
365 groupToGtype(i) = g
366 endif
367 enddo
368 endif
369 enddo
370
371 !! allocate the gtypeCutoffMap here.
372
373 !! then we do a double loop over all the group TYPES to find the cutoff
374 !! map between groups of two types
375
376 do i = 1, nGroupTypes
377 do j = 1, nGroupTypes
378
379 select case(cutoffPolicy)
380 case(TRADITIONAL_CUTOFF_POLICY)
381 thisRcut = maxval(gtypeMaxCutoff)
382 case(MIX_CUTOFF_POLICY)
383 thisRcut = 0.5_dp * (gtypeMaxCutoff(i) + gtypeMaxCutoff(j))
384 case(MAX_CUTOFF_POLICY)
385 thisRcut = max(gtypeMaxCutoff(i), gtypeMaxCutoff(j))
386 case default
387 call handleError("createGtypeCutoffMap", "Unknown Cutoff Policy")
388 return
389 end select
390 gtypeCutoffMap(i,j)%rcut = thisRcut
391 gtypeCutoffMap(i,j)%rcutsq = thisRcut*thisRcut
392 skin = defaultRlist - defaultRcut
393 gtypeCutoffMap(i,j)%rlistsq = (thisRcut + skin)**2
394 enddo
395 enddo
396
397 haveGtypeCutoffMap = .true.
398 end subroutine createGtypeCutoffMap
399
400 subroutine setDefaultCutoffs(defRcut, defRsw, defRlist, cutPolicy)
401 real(kind=dp),intent(in) :: defRcut, defRsw, defRlist
402 integer, intent(in) :: cutPolicy
403
404 defaultRcut = defRcut
405 defaultRsw = defRsw
406 defaultRlist = defRlist
407 cutoffPolicy = cutPolicy
408 end subroutine setDefaultCutoffs
409
410 subroutine setCutoffPolicy(cutPolicy)
411
412 integer, intent(in) :: cutPolicy
413 cutoffPolicy = cutPolicy
414 call createGtypeCutoffMap()
415
416 end subroutine setCutoffPolicy
417
418
419 subroutine setSimVariables()
420 SIM_uses_DirectionalAtoms = SimUsesDirectionalAtoms()
421 SIM_uses_EAM = SimUsesEAM()
422 SIM_uses_RF = SimUsesRF()
423 SIM_requires_postpair_calc = SimRequiresPostpairCalc()
424 SIM_requires_prepair_calc = SimRequiresPrepairCalc()
425 SIM_uses_PBC = SimUsesPBC()
426
427 haveSIMvariables = .true.
428
429 return
430 end subroutine setSimVariables
431
432 subroutine doReadyCheck(error)
433 integer, intent(out) :: error
434
435 integer :: myStatus
436
437 error = 0
438
439 if (.not. haveInteractionHash) then
440 myStatus = 0
441 call createInteractionHash(myStatus)
442 if (myStatus .ne. 0) then
443 write(default_error, *) 'createInteractionHash failed in doForces!'
444 error = -1
445 return
446 endif
447 endif
448
449 if (.not. haveGtypeCutoffMap) then
450 myStatus = 0
451 call createGtypeCutoffMap(myStatus)
452 if (myStatus .ne. 0) then
453 write(default_error, *) 'createGtypeCutoffMap failed in doForces!'
454 error = -1
455 return
456 endif
457 endif
458
459 if (.not. haveSIMvariables) then
460 call setSimVariables()
461 endif
462
463 if (.not. haveRlist) then
464 write(default_error, *) 'rList has not been set in doForces!'
465 error = -1
466 return
467 endif
468
469 if (.not. haveNeighborList) then
470 write(default_error, *) 'neighbor list has not been initialized in doForces!'
471 error = -1
472 return
473 end if
474
475 if (.not. haveSaneForceField) then
476 write(default_error, *) 'Force Field is not sane in doForces!'
477 error = -1
478 return
479 end if
480
481 #ifdef IS_MPI
482 if (.not. isMPISimSet()) then
483 write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
484 error = -1
485 return
486 endif
487 #endif
488 return
489 end subroutine doReadyCheck
490
491
492 subroutine init_FF(use_RF, use_UW, use_DW, thisStat)
493
494 logical, intent(in) :: use_RF
495 logical, intent(in) :: use_UW
496 logical, intent(in) :: use_DW
497 integer, intent(out) :: thisStat
498 integer :: my_status, nMatches
499 integer :: corrMethod
500 integer, pointer :: MatchList(:) => null()
501 real(kind=dp) :: rcut, rrf, rt, dielect
502
503 !! assume things are copacetic, unless they aren't
504 thisStat = 0
505
506 !! Fortran's version of a cast:
507 FF_uses_RF = use_RF
508
509 !! set the electrostatic correction method
510 if (use_UW) then
511 corrMethod = 1
512 elseif (use_DW) then
513 corrMethod = 2
514 else
515 corrMethod = 0
516 endif
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 FF_uses_RF setting makes sense
548
549 if (FF_uses_Dipoles) then
550 if (FF_uses_RF) then
551 dielect = getDielect()
552 call initialize_rf(dielect)
553 endif
554 else
555 if (FF_uses_RF) 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 (FF_uses_RF .and. SIM_uses_RF) 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 real ( kind = dp ) :: ebalance
1058 integer :: me_i, me_j
1059
1060 integer :: iHash
1061
1062 r = sqrt(rijsq)
1063 vpair = 0.0d0
1064 fpair(1:3) = 0.0d0
1065
1066 #ifdef IS_MPI
1067 me_i = atid_row(i)
1068 me_j = atid_col(j)
1069 #else
1070 me_i = atid(i)
1071 me_j = atid(j)
1072 #endif
1073
1074 iHash = InteractionHash(me_i, me_j)
1075
1076 if ( iand(iHash, LJ_PAIR).ne.0 ) then
1077 call do_lj_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, f, do_pot)
1078 endif
1079
1080 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1081 call doElectrostaticPair(i, j, d, r, rijsq, sw, vpair, fpair, &
1082 pot, eFrame, f, t, do_pot, corrMethod)
1083
1084 if (FF_uses_RF .and. SIM_uses_RF) then
1085
1086 ! CHECK ME (RF needs to know about all electrostatic types)
1087 call accumulate_rf(i, j, r, eFrame, sw)
1088 call rf_correct_forces(i, j, d, r, eFrame, sw, f, fpair)
1089 endif
1090
1091 endif
1092
1093 if ( iand(iHash, STICKY_PAIR).ne.0 ) then
1094 call do_sticky_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1095 pot, A, f, t, do_pot)
1096 endif
1097
1098 if ( iand(iHash, STICKYPOWER_PAIR).ne.0 ) then
1099 call do_sticky_power_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1100 pot, A, f, t, do_pot)
1101 endif
1102
1103 if ( iand(iHash, GAYBERNE_PAIR).ne.0 ) then
1104 call do_gb_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1105 pot, A, f, t, do_pot)
1106 endif
1107
1108 if ( iand(iHash, GAYBERNE_LJ).ne.0 ) then
1109 ! call do_gblj_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1110 ! pot, A, f, t, do_pot)
1111 endif
1112
1113 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1114 call do_eam_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, f, &
1115 do_pot)
1116 endif
1117
1118 if ( iand(iHash, SHAPE_PAIR).ne.0 ) then
1119 call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1120 pot, A, f, t, do_pot)
1121 endif
1122
1123 if ( iand(iHash, SHAPE_LJ).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 end subroutine do_pair
1129
1130 subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, &
1131 do_pot, do_stress, eFrame, A, f, t, pot)
1132
1133 real( kind = dp ) :: pot, sw
1134 real( kind = dp ), dimension(9,nLocal) :: eFrame
1135 real (kind=dp), dimension(9,nLocal) :: A
1136 real (kind=dp), dimension(3,nLocal) :: f
1137 real (kind=dp), dimension(3,nLocal) :: t
1138
1139 logical, intent(inout) :: do_pot, do_stress
1140 integer, intent(in) :: i, j
1141 real ( kind = dp ), intent(inout) :: rijsq, rcijsq
1142 real ( kind = dp ) :: r, rc
1143 real ( kind = dp ), intent(inout) :: d(3), dc(3)
1144
1145 integer :: me_i, me_j, iHash
1146
1147 #ifdef IS_MPI
1148 me_i = atid_row(i)
1149 me_j = atid_col(j)
1150 #else
1151 me_i = atid(i)
1152 me_j = atid(j)
1153 #endif
1154
1155 iHash = InteractionHash(me_i, me_j)
1156
1157 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1158 call calc_EAM_prepair_rho(i, j, d, r, rijsq )
1159 endif
1160
1161 end subroutine do_prepair
1162
1163
1164 subroutine do_preforce(nlocal,pot)
1165 integer :: nlocal
1166 real( kind = dp ) :: pot
1167
1168 if (FF_uses_EAM .and. SIM_uses_EAM) then
1169 call calc_EAM_preforce_Frho(nlocal,pot)
1170 endif
1171
1172
1173 end subroutine do_preforce
1174
1175
1176 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
1177
1178 real (kind = dp), dimension(3) :: q_i
1179 real (kind = dp), dimension(3) :: q_j
1180 real ( kind = dp ), intent(out) :: r_sq
1181 real( kind = dp ) :: d(3), scaled(3)
1182 integer i
1183
1184 d(1:3) = q_j(1:3) - q_i(1:3)
1185
1186 ! Wrap back into periodic box if necessary
1187 if ( SIM_uses_PBC ) then
1188
1189 if( .not.boxIsOrthorhombic ) then
1190 ! calc the scaled coordinates.
1191
1192 scaled = matmul(HmatInv, d)
1193
1194 ! wrap the scaled coordinates
1195
1196 scaled = scaled - anint(scaled)
1197
1198
1199 ! calc the wrapped real coordinates from the wrapped scaled
1200 ! coordinates
1201
1202 d = matmul(Hmat,scaled)
1203
1204 else
1205 ! calc the scaled coordinates.
1206
1207 do i = 1, 3
1208 scaled(i) = d(i) * HmatInv(i,i)
1209
1210 ! wrap the scaled coordinates
1211
1212 scaled(i) = scaled(i) - anint(scaled(i))
1213
1214 ! calc the wrapped real coordinates from the wrapped scaled
1215 ! coordinates
1216
1217 d(i) = scaled(i)*Hmat(i,i)
1218 enddo
1219 endif
1220
1221 endif
1222
1223 r_sq = dot_product(d,d)
1224
1225 end subroutine get_interatomic_vector
1226
1227 subroutine zero_work_arrays()
1228
1229 #ifdef IS_MPI
1230
1231 q_Row = 0.0_dp
1232 q_Col = 0.0_dp
1233
1234 q_group_Row = 0.0_dp
1235 q_group_Col = 0.0_dp
1236
1237 eFrame_Row = 0.0_dp
1238 eFrame_Col = 0.0_dp
1239
1240 A_Row = 0.0_dp
1241 A_Col = 0.0_dp
1242
1243 f_Row = 0.0_dp
1244 f_Col = 0.0_dp
1245 f_Temp = 0.0_dp
1246
1247 t_Row = 0.0_dp
1248 t_Col = 0.0_dp
1249 t_Temp = 0.0_dp
1250
1251 pot_Row = 0.0_dp
1252 pot_Col = 0.0_dp
1253 pot_Temp = 0.0_dp
1254
1255 rf_Row = 0.0_dp
1256 rf_Col = 0.0_dp
1257 rf_Temp = 0.0_dp
1258
1259 #endif
1260
1261 if (FF_uses_EAM .and. SIM_uses_EAM) then
1262 call clean_EAM()
1263 endif
1264
1265 rf = 0.0_dp
1266 tau_Temp = 0.0_dp
1267 virial_Temp = 0.0_dp
1268 end subroutine zero_work_arrays
1269
1270 function skipThisPair(atom1, atom2) result(skip_it)
1271 integer, intent(in) :: atom1
1272 integer, intent(in), optional :: atom2
1273 logical :: skip_it
1274 integer :: unique_id_1, unique_id_2
1275 integer :: me_i,me_j
1276 integer :: i
1277
1278 skip_it = .false.
1279
1280 !! there are a number of reasons to skip a pair or a particle
1281 !! mostly we do this to exclude atoms who are involved in short
1282 !! range interactions (bonds, bends, torsions), but we also need
1283 !! to exclude some overcounted interactions that result from
1284 !! the parallel decomposition
1285
1286 #ifdef IS_MPI
1287 !! in MPI, we have to look up the unique IDs for each atom
1288 unique_id_1 = AtomRowToGlobal(atom1)
1289 #else
1290 !! in the normal loop, the atom numbers are unique
1291 unique_id_1 = atom1
1292 #endif
1293
1294 !! We were called with only one atom, so just check the global exclude
1295 !! list for this atom
1296 if (.not. present(atom2)) then
1297 do i = 1, nExcludes_global
1298 if (excludesGlobal(i) == unique_id_1) then
1299 skip_it = .true.
1300 return
1301 end if
1302 end do
1303 return
1304 end if
1305
1306 #ifdef IS_MPI
1307 unique_id_2 = AtomColToGlobal(atom2)
1308 #else
1309 unique_id_2 = atom2
1310 #endif
1311
1312 #ifdef IS_MPI
1313 !! this situation should only arise in MPI simulations
1314 if (unique_id_1 == unique_id_2) then
1315 skip_it = .true.
1316 return
1317 end if
1318
1319 !! this prevents us from doing the pair on multiple processors
1320 if (unique_id_1 < unique_id_2) then
1321 if (mod(unique_id_1 + unique_id_2,2) == 0) then
1322 skip_it = .true.
1323 return
1324 endif
1325 else
1326 if (mod(unique_id_1 + unique_id_2,2) == 1) then
1327 skip_it = .true.
1328 return
1329 endif
1330 endif
1331 #endif
1332
1333 !! the rest of these situations can happen in all simulations:
1334 do i = 1, nExcludes_global
1335 if ((excludesGlobal(i) == unique_id_1) .or. &
1336 (excludesGlobal(i) == unique_id_2)) then
1337 skip_it = .true.
1338 return
1339 endif
1340 enddo
1341
1342 do i = 1, nSkipsForAtom(atom1)
1343 if (skipsForAtom(atom1, i) .eq. unique_id_2) then
1344 skip_it = .true.
1345 return
1346 endif
1347 end do
1348
1349 return
1350 end function skipThisPair
1351
1352 function FF_UsesDirectionalAtoms() result(doesit)
1353 logical :: doesit
1354 doesit = FF_uses_DirectionalAtoms
1355 end function FF_UsesDirectionalAtoms
1356
1357 function FF_RequiresPrepairCalc() result(doesit)
1358 logical :: doesit
1359 doesit = FF_uses_EAM
1360 end function FF_RequiresPrepairCalc
1361
1362 function FF_RequiresPostpairCalc() result(doesit)
1363 logical :: doesit
1364 doesit = FF_uses_RF
1365 end function FF_RequiresPostpairCalc
1366
1367 #ifdef PROFILE
1368 function getforcetime() result(totalforcetime)
1369 real(kind=dp) :: totalforcetime
1370 totalforcetime = forcetime
1371 end function getforcetime
1372 #endif
1373
1374 !! This cleans componets of force arrays belonging only to fortran
1375
1376 subroutine add_stress_tensor(dpair, fpair)
1377
1378 real( kind = dp ), dimension(3), intent(in) :: dpair, fpair
1379
1380 ! because the d vector is the rj - ri vector, and
1381 ! because fx, fy, fz are the force on atom i, we need a
1382 ! negative sign here:
1383
1384 tau_Temp(1) = tau_Temp(1) - dpair(1) * fpair(1)
1385 tau_Temp(2) = tau_Temp(2) - dpair(1) * fpair(2)
1386 tau_Temp(3) = tau_Temp(3) - dpair(1) * fpair(3)
1387 tau_Temp(4) = tau_Temp(4) - dpair(2) * fpair(1)
1388 tau_Temp(5) = tau_Temp(5) - dpair(2) * fpair(2)
1389 tau_Temp(6) = tau_Temp(6) - dpair(2) * fpair(3)
1390 tau_Temp(7) = tau_Temp(7) - dpair(3) * fpair(1)
1391 tau_Temp(8) = tau_Temp(8) - dpair(3) * fpair(2)
1392 tau_Temp(9) = tau_Temp(9) - dpair(3) * fpair(3)
1393
1394 virial_Temp = virial_Temp + &
1395 (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
1396
1397 end subroutine add_stress_tensor
1398
1399 end module doForces