ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/doForces.F90
Revision: 571
Committed: Tue Aug 9 22:33:37 2005 UTC (19 years, 9 months ago) by gezelter
File size: 36281 byte(s)
Log Message:
Breaky Breaky BREAKY breaky breaky

File Contents

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