ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 1214
Committed: Tue Jun 1 18:42:58 2004 UTC (20 years, 3 months ago) by gezelter
File size: 34726 byte(s)
Log Message:
Cutoff Groups for MPI

File Contents

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