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