15 double precision,
public ::
vc_mu = 1.d0
18 integer,
private,
parameter :: rho_ = 1
21 integer,
allocatable,
private,
protected :: mom(:)
24 integer,
private,
protected :: e_
44 character(len=*),
intent(in) :: files(:)
50 open(
unitpar, file=trim(files(n)), status=
"old")
60 integer,
intent(inout) :: phys_wider_stencil
61 logical,
intent(inout) :: phys_req_diagonal
82 phys_wider_stencil = 1
83 phys_req_diagonal = .true.
89 energy,qsourcesplit,active)
96 integer,
intent(in) :: ixI^L, ixO^L
97 double precision,
intent(in) :: qdt, x(ixI^S,1:ndim), wCT(ixI^S,1:nw)
98 double precision,
intent(inout) :: w(ixI^S,1:nw)
99 logical,
intent(in) :: energy,qsourcesplit
100 logical,
intent(inout) :: active
102 double precision:: lambda(ixI^S,ndir,ndir),tmp(ixI^S),tmp2(ixI^S),v(ixI^S,ndir),vlambda(ixI^S,ndir)
103 integer:: ix^L,idim,idir,jdir,iw
107 if(qsourcesplit .eqv.
vc_split)
then
114 if({ iximin^
d>ixmin^
d .or. iximax^
d<ixmax^
d|.or.})&
115 call mpistop(
"error for viscous source addition, 2 layers needed")
120 if({ iximin^
d>ixmin^
d .or. iximax^
d<ixmax^
d|.or.})&
121 call mpistop(
"error for viscous source addition"//&
122 "requested fourth order gradients: 4 layers needed")
128 v(ixi^s,idir)=wct(ixi^s,mom(idir))/wct(ixi^s,rho_)
136 do idim=1,ndim;
do idir=1,ndir
139 tmp(ixi^s)=v(ixi^s,idir)
143 if (idim==
r_ .and. (idir==2 .or. idir==
phi_))
then
144 tmp(ixi^s) = tmp(ixi^s)/x(ixi^s,1)
146 elseif (idim==2 .and. idir==
phi_)
then
147 tmp(ixi^s)=tmp(ixi^s)/dsin(x(ixi^s,2))
151 call gradient(tmp,ixi^l,ix^l,idim,tmp2)
156 if (idim==
r_ .and. (idir==2 .or. idir==
phi_))
then
157 tmp2(ix^s) = tmp2(ix^s)*x(ix^s,1)
159 elseif (idim==2 .and. idir==
phi_ )
then
160 tmp2(ix^s)=tmp2(ix^s)*dsin(x(ix^s,2))
161 elseif (idim==2 .and. idir==2 )
then
162 tmp2(ix^s)=tmp2(ix^s)+v(ix^s,
r_)/x(ix^s,1)
163 elseif (idim==
phi_.and. idir==
phi_)
then
164 tmp2(ix^s)=tmp2(ix^s)+v(ix^s,
r_)/x(ix^s,1)+v(ix^s,2)/(x(ix^s,1)*dtan(x(ix^s,2)))
168 lambda(ix^s,idim,idir)= lambda(ix^s,idim,idir)+ tmp2(ix^s)
169 lambda(ix^s,idir,idim)= lambda(ix^s,idir,idim)+ tmp2(ix^s)
173 lambda(ix^s,1:ndir,1:ndir)=lambda(ix^s,1:ndir,1:ndir)*
vc_mu*qdt
181 tmp(ix^s)=tmp(ix^s)+lambda(ix^s,idir,idir)
183 tmp(ix^s)=tmp(ix^s)/3.d0
187 lambda(ix^s,idir,idir)=lambda(ix^s,idir,idir)-tmp(ix^s)
194 tmp(ix^s)=lambda(ix^s,idir,idim)
199 if (idim==
r_ .and. idir==
r_ ) tmp(ix^s) = tmp(ix^s)*x(ix^s,1)**two
200 if (idim==
r_ .and. (idir==2 .or. idir==
phi_)) tmp(ix^s) = tmp(ix^s)*x(ix^s,1)**3.d0
202 if (idim==2 .and. (idir==
r_ .or. idir==2)) tmp(ix^s) = tmp(ix^s)*dsin(x(ix^s,2))
203 if (idim==2 .and. idir==
phi_ ) tmp(ix^s) = tmp(ix^s)*dsin(x(ix^s,2))**two
206 call gradient(tmp,ixi^l,ixo^l,idim,tmp2)
211 if (idim==
r_ .and. idir==
r_ ) tmp2(ixo^s) = tmp2(ixo^s)/(x(ixo^s,1)**two)
212 if (idim==
r_ .and. (idir==2 .or. idir==
phi_)) tmp2(ixo^s) = tmp2(ixo^s)/(x(ixo^s,1)**3.d0)
214 if (idim==2 .and. (idir==
r_ .or. idir==2)) tmp2(ixo^s) = tmp2(ixo^s)/(dsin(x(ixo^s,2)))
215 if (idim==2 .and. idir==
phi_ ) tmp2(ixo^s) = tmp2(ixo^s)/(dsin(x(ixo^s,2))**two)
218 w(ixo^s,mom(idir))=w(ixo^s,mom(idir))+tmp2(ixo^s)
224 if (
coordinate==
spherical .and. idir==2 ) w(ixo^s,mom(idir))=w(ixo^s,mom(idir))-lambda(ixo^s,
phi_,
phi_)/(x(ixo^s,1)/dtan(x(ixo^s,2)))
234 vlambda(ixi^s,idim)=vlambda(ixi^s,idim)+v(ixi^s,idir)*lambda(ixi^s,idir,idim)
238 w(ixo^s,e_)=w(ixo^s,e_)+tmp2(ixo^s)
248 integer,
intent(in) :: ixI^L, ixO^L
249 double precision,
intent(in) :: dx^D, x(ixI^S,1:ndim)
250 double precision,
intent(in) :: w(ixI^S,1:nw)
251 double precision,
intent(inout) :: dtnew
253 double precision :: tmp(ixI^S)
254 double precision:: dtdiff_visc, dxinv2(1:ndim)
259 tmp(ixo^s)=
vc_mu/w(ixo^s,rho_)
261 ^d&dxinv2(^d)=one/dx^d**2;
263 dtdiff_visc=
dtdiffpar/maxval(tmp(ixo^s)*dxinv2(idim))
265 dtnew=min(dtnew,dtdiff_visc)
279 integer,
intent(in) :: ixi^
l, ixo^
l, idim
280 double precision,
intent(in) :: w(ixi^s, 1:nw), x(ixi^s, 1:^nd)
281 double precision,
intent(inout) :: f(ixi^s, nwflux)
282 logical,
intent(in) :: energy
284 double precision :: v(ixi^s,1:
ndir)
286 double precision :: divergence(ixi^s)
288 double precision:: lambda(ixi^s,
ndir)
293 v(ixi^s,i)=w(ixi^s,i+1)
298 lambda(ixo^s,idim) = lambda(ixo^s,idim) - (2.d0/3.d0) * divergence(ixo^s)
302 f(ixo^s, mom(idir)) = f(ixo^s, mom(idir)) -
vc_mu * lambda(ixo^s,idir)
303 if (energy) f(ixo^s, e_) = f(ixo^s, e_) -
vc_mu * lambda(ixo^s,idir) * v(ixi^s,idir)
313 integer,
intent(in) :: ixI^L, ixO^L, idim
314 double precision,
intent(in) :: w(ixI^S, 1:nw), x(ixI^S, 1:ndim)
315 double precision,
intent(out) :: cross(ixI^S,ndir)
317 double precision :: tmp(ixI^S), v(ixI^S)
320 if (ndir/=ndim)
call mpistop(
"This formula are probably wrong for ndim/=ndir")
334 v(ixi^s)=w(ixi^s,mom(1))
336 cross(ixi^s,2)=tmp(ixi^s)
337 v(ixi^s)=w(ixi^s,mom(2))/x(ixi^s,1)
339 cross(ixi^s,2)=cross(ixi^s,2)+tmp(ixi^s)*x(ixi^s,1)
341 elseif (idim==2)
then
343 v(ixi^s)=w(ixi^s,mom(1))
346 cross(ixi^s,1)=tmp(ixi^s)
347 v(ixi^s)=w(ixi^s,mom(2))/x(ixi^s,1)
349 cross(ixi^s,1)=cross(ixi^s,1)+tmp(ixi^s)*x(ixi^s,1)
351 v(ixi^s)=w(ixi^s,mom(2))
353 cross(ixi^s,2)=two*tmp(ixi^s)
354 v(ixi^s)=w(ixi^s,mom(1))
355 cross(ixi^s,2)=cross(ixi^s,2)+two*v(ixi^s)/x(ixi^s,1)
357 v(ixi^s)=w(ixi^s,mom(3))
360 cross(ixi^s,3)=tmp(ixi^s)
361 v(ixi^s)=w(ixi^s,mom(2))
365 cross(ixi^s,3)=cross(ixi^s,3)+tmp(ixi^s)
367 elseif (idim==3)
then
372 v(ixi^s)=w(ixi^s,mom(3))
374 cross(ixi^s,2)=tmp(ixi^s)
375 v(ixi^s)=w(ixi^s,mom(2))
377 cross(ixi^s,2)=cross(ixi^s,2)+tmp(ixi^s)
383 v(ixi^s)=w(ixi^s,mom(1))
385 cross(ixi^s,1)=two*tmp(ixi^s)
387 v(ixi^s)=w(ixi^s,mom(1))
390 cross(ixi^s,2)=tmp(ixi^s)
391 v(ixi^s)=w(ixi^s,mom(2))/x(ixi^s,1)
393 cross(ixi^s,2)=cross(ixi^s,2)+tmp(ixi^s)*x(ixi^s,1)
397 v(ixi^s)=w(ixi^s,mom(1))
399 cross(ixi^s,3)=tmp(ixi^s)
400 v(ixi^s)=w(ixi^s,mom(3))/x(ixi^s,1)
402 cross(ixi^s,3)=cross(ixi^s,3)+tmp(ixi^s)*x(ixi^s,1)
404 elseif (idim==2)
then
406 v(ixi^s)=w(ixi^s,mom(1))
409 cross(ixi^s,1)=tmp(ixi^s)
410 v(ixi^s)=w(ixi^s,mom(2))/x(ixi^s,1)
412 cross(ixi^s,1)=cross(ixi^s,1)+tmp(ixi^s)*x(ixi^s,1)
414 v(ixi^s)=w(ixi^s,mom(2))
416 cross(ixi^s,2)=two*tmp(ixi^s)
417 v(ixi^s)=w(ixi^s,mom(1))
418 cross(ixi^s,2)=cross(ixi^s,2)+two*v(ixi^s)/x(ixi^s,1)
422 v(ixi^s)=w(ixi^s,mom(2))
424 cross(ixi^s,3)=tmp(ixi^s)
425 v(ixi^s)=w(ixi^s,mom(3))/dsin(x(ixi^s,2))
427 cross(ixi^s,3)=cross(ixi^s,3)+tmp(ixi^s)*dsin(x(ixi^s,2))
430 elseif (idim==3)
then
432 v(ixi^s)=w(ixi^s,mom(1))
434 cross(ixi^s,1)=tmp(ixi^s)
435 v(ixi^s)=w(ixi^s,mom(3))/x(ixi^s,1)
437 cross(ixi^s,1)=cross(ixi^s,1)+tmp(ixi^s)*x(ixi^s,1)
439 v(ixi^s)=w(ixi^s,mom(2))
441 cross(ixi^s,2)=tmp(ixi^s)
442 v(ixi^s)=w(ixi^s,mom(3))/dsin(x(ixi^s,2))
444 cross(ixi^s,2)=cross(ixi^s,2)+tmp(ixi^s)*dsin(x(ixi^s,2))
446 v(ixi^s)=w(ixi^s,mom(3))
448 cross(ixi^s,3)=two*tmp(ixi^s)
449 v(ixi^s)=w(ixi^s,mom(1))
450 cross(ixi^s,3)=cross(ixi^s,3)+two*v(ixi^s)/x(ixi^s,1)
451 v(ixi^s)=w(ixi^s,mom(2))
452 cross(ixi^s,3)=cross(ixi^s,3)+two*v(ixi^s)/(x(ixi^s,1)*dtan(x(ixi^s,2)))
456 call mpistop(
"Unknown geometry specified")
466 integer,
intent(in) :: ixI^L, ixO^L, idim
467 double precision,
intent(in) :: w(ixI^S, 1:nw), x(ixI^S, 1:^ND)
468 double precision,
intent(out) :: cross(ixI^S,ndir)
470 double precision :: tmp(ixI^S), v(ixI^S)
473 v(ixi^s)=w(ixi^s,mom(idim))
475 call gradient(v,ixi^l,ixo^l,idir,tmp)
476 cross(ixo^s,idir)=tmp(ixo^s)
479 v(ixi^s)=w(ixi^s,mom(idir))
480 call gradient(v,ixi^l,ixo^l,idim,tmp)
481 cross(ixo^s,idir)=cross(ixo^s,idir)+tmp(ixo^s)
490 integer,
intent(in) :: ixI^L, ixO^L
491 double precision,
intent(in) :: qdt, x(ixI^S, 1:ndim)
492 double precision,
intent(inout) :: wCT(ixI^S, 1:nw), w(ixI^S, 1:nw)
496 double precision :: vv(ixI^S), divergence(ixI^S)
497 double precision :: tmp(ixI^S),tmp1(ixI^S)
507 call gradient(wct(ixi^s,mom(2)),ixi^l,ixo^l,2,tmp1)
508 tmp(ixo^s)=two*(tmp1(ixo^s)+wct(ixo^s,mom(1))/x(ixo^s,1))
510 call divvector(wct(ixi^s,mom(1:
ndir)),ixi^l,ixo^l,divergence)
511 tmp(ixo^s) = tmp(ixo^s) - (2.d0/3.d0) * divergence(ixo^s)
513 w(ixo^s,mom(1))=w(ixo^s,mom(1))-qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
515 call gradient(wct(ixi^s,mom(1)),ixi^l,ixo^l,2,tmp)
516 vv(ixi^s)=wct(ixi^s,mom(2))/x(ixi^s,1)
517 call gradient(vv,ixi^l,ixo^l,1,tmp1)
518 tmp(ixo^s)=tmp(ixo^s)+tmp1(ixo^s)*x(ixo^s,1)
520 w(ixo^s,mom(2))=w(ixo^s,mom(2))+qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
526 call gradient(wct(ixi^s,mom(2)),ixi^l,ixo^l,2,tmp1)
527 tmp(ixo^s)=two*(tmp1(ixo^s)+wct(ixo^s,mom(1))/x(ixo^s,1))
529 call divvector(wct(ixi^s,mom(1:
ndir)),ixi^l,ixo^l,divergence)
530 tmp(ixo^s) = tmp(ixo^s) - (2.d0/3.d0) * divergence(ixo^s)
532 w(ixo^s,mom(1))=w(ixo^s,mom(1))-qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
537 call gradient(wct(ixi^s,mom(3)),ixi^l,ixo^l,3,tmp1)
538 tmp(ixo^s)=two*(tmp1(ixo^s)+wct(ixo^s,mom(1))/x(ixo^s,1)+wct(ixo^s,mom(2))/(x(ixo^s,1)*dtan(x(ixo^s,2))))
541 tmp(ixo^s) = tmp(ixo^s) - (2.d0/3.d0) * divergence(ixo^s)
543 w(ixo^s,mom(1))=w(ixo^s,mom(1))-qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
546 w(ixo^s,mom(2))=w(ixo^s,mom(2))-qdt*
vc_mu*tmp(ixo^s)/(x(ixo^s,1)*dtan(x(ixo^s,2)))
549 call gradient(wct(ixi^s,mom(1)),ixi^l,ixo^l,2,tmp)
550 vv(ixi^s)=wct(ixi^s,mom(2))/x(ixi^s,1)
551 call gradient(vv,ixi^l,ixo^l,1,tmp1)
552 tmp(ixo^s)=tmp(ixo^s)+tmp1(ixo^s)*x(ixo^s,1)
554 w(ixo^s,mom(2))=w(ixo^s,mom(2))+qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
557 call gradient(wct(ixi^s,mom(1)),ixi^l,ixo^l,3,tmp)
559 vv(ixi^s)=wct(ixi^s,mom(3))/x(ixi^s,1)
560 call gradient(vv,ixi^l,ixo^l,1,tmp1)
561 tmp(ixo^s)=tmp(ixo^s)+tmp1(ixo^s)*x(ixo^s,1)
563 w(ixo^s,mom(3))=w(ixo^s,mom(3))+qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
566 call gradient(wct(ixi^s,mom(2)),ixi^l,ixo^l,3,tmp)
569 vv(ixi^s)=wct(ixi^s,mom(3))/dsin(x(ixi^s,2))
570 call gradient(vv,ixi^l,ixo^l,2,tmp1)
571 tmp(ixo^s)=tmp(ixo^s)+tmp1(ixo^s)*dsin(x(ixo^s,2))
573 w(ixo^s,mom(3))=w(ixo^s,mom(3))+qdt*
vc_mu*tmp(ixo^s)/(x(ixo^s,1)*dtan(x(ixo^s,2)))
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module with geometry-related routines (e.g., divergence, curl)
integer, parameter spherical
integer, parameter cartesian
integer, parameter cylindrical
integer, parameter cartesian_stretched
subroutine gradient(q, ixIL, ixOL, idir, gradq)
Calculate gradient of a scalar q within ixL in direction idir.
subroutine divvector(qvec, ixIL, ixOL, divq, fourthorder, sixthorder)
Calculate divergence of a vector qvec within ixL.
This module contains definitions of global parameters and variables and some generic functions/subrou...
double precision dtdiffpar
For resistive MHD, the time step is also limited by the diffusion time: .
integer, parameter unitpar
file handle for IO
logical any_source_split
if any normal source term is added in split fasion
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
double precision, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
integer r_
Indices for cylindrical coordinates FOR TESTS, negative value when not used:
The module add viscous source terms and check time step.
subroutine get_crossgrad(ixIL, ixOL, x, w, idim, cross)
subroutine viscosity_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
subroutine viscosity_add_source(qdt, ixIL, ixOL, wCT, w, x, energy, qsourcesplit, active)
subroutine viscosity_init(phys_wider_stencil, phys_req_diagonal)
Initialize the module.
subroutine cart_cross_grad(ixIL, ixOL, x, w, idim, cross)
yields d_i v_j + d_j v_i for a given i, OK in Cartesian and for some tensor terms in cylindrical (rr ...
subroutine visc_add_source_geom(qdt, ixIL, ixOL, wCT, w, x)
double precision, public vc_mu
Viscosity coefficient.
logical vc_4th_order
fourth order
logical vc_split
source split or not
subroutine, public visc_get_flux_prim(w, x, ixIL, ixOL, idim, f, energy)
subroutine vc_params_read(files)
Read this module"s parameters from a file.
logical viscindiv
whether to compute the viscous terms as fluxes (ie in the div on the LHS), or not (by default)