ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2280
Committed: Thu Sep 1 20:17:55 2005 UTC (18 years, 10 months ago) by gezelter
File size: 41725 byte(s)
Log Message:
wrote createGtypeCutoffMap

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