15 double precision,
public ::
vc_mu = 1.d0
28 integer,
private,
parameter :: rho_ = 1
31 integer,
allocatable,
private,
protected :: mom(:)
34 integer,
private,
protected :: e_
43 character(len=*),
intent(in) :: files(:)
49 open(
unitpar, file=trim(files(n)), status=
"old")
59 integer,
intent(inout) :: phys_wider_stencil
60 logical,
intent(inout) :: phys_req_diagonal
81 phys_wider_stencil = 1
82 phys_req_diagonal = .true.
88 energy,qsourcesplit,active)
95 integer,
intent(in) :: ixI^L, ixO^L
96 double precision,
intent(in) :: qdt, x(ixI^S,1:ndim), wCT(ixI^S,1:nw)
97 double precision,
intent(inout) :: w(ixI^S,1:nw)
98 logical,
intent(in) :: energy,qsourcesplit
99 logical,
intent(inout) :: active
101 integer:: ix^L,idim,idir,jdir,iw
102 double precision:: lambda(ixI^S,ndir,ndir),tmp(ixI^S),tmp2(ixI^S),v(ixI^S,ndir),vlambda(ixI^S,ndir)
106 if(qsourcesplit .eqv.
vc_split)
then
113 if({ iximin^
d>ixmin^
d .or. iximax^
d<ixmax^
d|.or.})&
114 call mpistop(
"error for viscous source addition, 2 layers needed")
119 if({ iximin^
d>ixmin^
d .or. iximax^
d<ixmax^
d|.or.})&
120 call mpistop(
"error for viscous source addition"//&
121 "requested fourth order gradients: 4 layers needed")
127 v(ixi^s,idir)=wct(ixi^s,mom(idir))/wct(ixi^s,rho_)
135 do idim=1,ndim;
do idir=1,ndir
138 tmp(ixi^s)=v(ixi^s,idir)
142 if (idim==
r_ .and. (idir==2 .or. idir==
phi_))
then
143 tmp(ixi^s) = tmp(ixi^s)/x(ixi^s,1)
145 elseif (idim==2 .and. idir==
phi_)
then
146 tmp(ixi^s)=tmp(ixi^s)/dsin(x(ixi^s,2))
150 call gradient(tmp,ixi^l,ix^l,idim,tmp2)
155 if (idim==
r_ .and. (idir==2 .or. idir==
phi_))
then
156 tmp2(ix^s) = tmp2(ix^s)*x(ix^s,1)
158 elseif (idim==2 .and. idir==
phi_ )
then
159 tmp2(ix^s)=tmp2(ix^s)*dsin(x(ix^s,2))
160 elseif (idim==2 .and. idir==2 )
then
161 tmp2(ix^s)=tmp2(ix^s)+v(ix^s,
r_)/x(ix^s,1)
162 elseif (idim==
phi_.and. idir==
phi_)
then
163 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)))
167 lambda(ix^s,idim,idir)= lambda(ix^s,idim,idir)+ tmp2(ix^s)
168 lambda(ix^s,idir,idim)= lambda(ix^s,idir,idim)+ tmp2(ix^s)
172 lambda(ix^s,1:ndir,1:ndir)=lambda(ix^s,1:ndir,1:ndir)*
vc_mu*qdt
180 tmp(ix^s)=tmp(ix^s)+lambda(ix^s,idir,idir)
182 tmp(ix^s)=tmp(ix^s)/3.d0
186 lambda(ix^s,idir,idir)=lambda(ix^s,idir,idir)-tmp(ix^s)
193 tmp(ix^s)=lambda(ix^s,idir,idim)
198 if (idim==
r_ .and. idir==
r_ ) tmp(ix^s) = tmp(ix^s)*x(ix^s,1)**two
199 if (idim==
r_ .and. (idir==2 .or. idir==
phi_)) tmp(ix^s) = tmp(ix^s)*x(ix^s,1)**3.d0
201 if (idim==2 .and. (idir==
r_ .or. idir==2)) tmp(ix^s) = tmp(ix^s)*dsin(x(ix^s,2))
202 if (idim==2 .and. idir==
phi_ ) tmp(ix^s) = tmp(ix^s)*dsin(x(ix^s,2))**two
205 call gradient(tmp,ixi^l,ixo^l,idim,tmp2)
210 if (idim==
r_ .and. idir==
r_ ) tmp2(ixo^s) = tmp2(ixo^s)/(x(ixo^s,1)**two)
211 if (idim==
r_ .and. (idir==2 .or. idir==
phi_)) tmp2(ixo^s) = tmp2(ixo^s)/(x(ixo^s,1)**3.d0)
213 if (idim==2 .and. (idir==
r_ .or. idir==2)) tmp2(ixo^s) = tmp2(ixo^s)/(dsin(x(ixo^s,2)))
214 if (idim==2 .and. idir==
phi_ ) tmp2(ixo^s) = tmp2(ixo^s)/(dsin(x(ixo^s,2))**two)
217 w(ixo^s,mom(idir))=w(ixo^s,mom(idir))+tmp2(ixo^s)
223 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)))
233 vlambda(ixi^s,idim)=vlambda(ixi^s,idim)+v(ixi^s,idir)*lambda(ixi^s,idir,idim)
237 w(ixo^s,e_)=w(ixo^s,e_)+tmp2(ixo^s)
247 integer,
intent(in) :: ixI^L, ixO^L
248 double precision,
intent(in) :: dx^D, x(ixI^S,1:ndim)
249 double precision,
intent(in) :: w(ixI^S,1:nw)
250 double precision,
intent(inout) :: dtnew
252 double precision :: tmp(ixI^S)
253 double precision:: dtdiff_visc, dxinv2(1:ndim)
258 tmp(ixo^s)=
vc_mu/w(ixo^s,rho_)
260 ^d&dxinv2(^d)=one/dx^d**2;
262 dtdiff_visc=
dtdiffpar/maxval(tmp(ixo^s)*dxinv2(idim))
264 dtnew=min(dtnew,dtdiff_visc)
278 integer,
intent(in) :: ixi^
l, ixo^
l, idim
279 double precision,
intent(in) :: w(ixi^s, 1:nw), x(ixi^s, 1:^nd)
280 double precision,
intent(inout) :: f(ixi^s, nwflux)
281 logical,
intent(in) :: energy
283 double precision :: v(ixi^s,1:
ndir)
285 double precision :: divergence(ixi^s)
287 double precision:: lambda(ixi^s,
ndir)
292 v(ixi^s,i)=w(ixi^s,i+1)
297 lambda(ixo^s,idim) = lambda(ixo^s,idim) - (2.d0/3.d0) * divergence(ixo^s)
301 f(ixo^s, mom(idir)) = f(ixo^s, mom(idir)) -
vc_mu * lambda(ixo^s,idir)
302 if (energy) f(ixo^s, e_) = f(ixo^s, e_) -
vc_mu * lambda(ixo^s,idir) * v(ixi^s,idir)
312 integer,
intent(in) :: ixI^L, ixO^L, idim
313 double precision,
intent(in) :: w(ixI^S, 1:nw), x(ixI^S, 1:ndim)
314 double precision,
intent(out) :: cross(ixI^S,ndir)
316 double precision :: tmp(ixI^S), v(ixI^S)
318 if (ndir/=ndim)
call mpistop(
"This formula are probably wrong for ndim/=ndir")
332 v(ixi^s)=w(ixi^s,mom(1))
334 cross(ixi^s,2)=tmp(ixi^s)
335 v(ixi^s)=w(ixi^s,mom(2))/x(ixi^s,1)
337 cross(ixi^s,2)=cross(ixi^s,2)+tmp(ixi^s)*x(ixi^s,1)
339 elseif (idim==2)
then
341 v(ixi^s)=w(ixi^s,mom(1))
344 cross(ixi^s,1)=tmp(ixi^s)
345 v(ixi^s)=w(ixi^s,mom(2))/x(ixi^s,1)
347 cross(ixi^s,1)=cross(ixi^s,1)+tmp(ixi^s)*x(ixi^s,1)
349 v(ixi^s)=w(ixi^s,mom(2))
351 cross(ixi^s,2)=two*tmp(ixi^s)
352 v(ixi^s)=w(ixi^s,mom(1))
353 cross(ixi^s,2)=cross(ixi^s,2)+two*v(ixi^s)/x(ixi^s,1)
355 v(ixi^s)=w(ixi^s,mom(3))
358 cross(ixi^s,3)=tmp(ixi^s)
359 v(ixi^s)=w(ixi^s,mom(2))
363 cross(ixi^s,3)=cross(ixi^s,3)+tmp(ixi^s)
365 elseif (idim==3)
then
370 v(ixi^s)=w(ixi^s,mom(3))
372 cross(ixi^s,2)=tmp(ixi^s)
373 v(ixi^s)=w(ixi^s,mom(2))
375 cross(ixi^s,2)=cross(ixi^s,2)+tmp(ixi^s)
381 v(ixi^s)=w(ixi^s,mom(1))
383 cross(ixi^s,1)=two*tmp(ixi^s)
385 v(ixi^s)=w(ixi^s,mom(1))
388 cross(ixi^s,2)=tmp(ixi^s)
389 v(ixi^s)=w(ixi^s,mom(2))/x(ixi^s,1)
391 cross(ixi^s,2)=cross(ixi^s,2)+tmp(ixi^s)*x(ixi^s,1)
395 v(ixi^s)=w(ixi^s,mom(1))
397 cross(ixi^s,3)=tmp(ixi^s)
398 v(ixi^s)=w(ixi^s,mom(3))/x(ixi^s,1)
400 cross(ixi^s,3)=cross(ixi^s,3)+tmp(ixi^s)*x(ixi^s,1)
402 elseif (idim==2)
then
404 v(ixi^s)=w(ixi^s,mom(1))
407 cross(ixi^s,1)=tmp(ixi^s)
408 v(ixi^s)=w(ixi^s,mom(2))/x(ixi^s,1)
410 cross(ixi^s,1)=cross(ixi^s,1)+tmp(ixi^s)*x(ixi^s,1)
412 v(ixi^s)=w(ixi^s,mom(2))
414 cross(ixi^s,2)=two*tmp(ixi^s)
415 v(ixi^s)=w(ixi^s,mom(1))
416 cross(ixi^s,2)=cross(ixi^s,2)+two*v(ixi^s)/x(ixi^s,1)
420 v(ixi^s)=w(ixi^s,mom(2))
422 cross(ixi^s,3)=tmp(ixi^s)
423 v(ixi^s)=w(ixi^s,mom(3))/dsin(x(ixi^s,2))
425 cross(ixi^s,3)=cross(ixi^s,3)+tmp(ixi^s)*dsin(x(ixi^s,2))
428 elseif (idim==3)
then
430 v(ixi^s)=w(ixi^s,mom(1))
432 cross(ixi^s,1)=tmp(ixi^s)
433 v(ixi^s)=w(ixi^s,mom(3))/x(ixi^s,1)
435 cross(ixi^s,1)=cross(ixi^s,1)+tmp(ixi^s)*x(ixi^s,1)
437 v(ixi^s)=w(ixi^s,mom(2))
439 cross(ixi^s,2)=tmp(ixi^s)
440 v(ixi^s)=w(ixi^s,mom(3))/dsin(x(ixi^s,2))
442 cross(ixi^s,2)=cross(ixi^s,2)+tmp(ixi^s)*dsin(x(ixi^s,2))
444 v(ixi^s)=w(ixi^s,mom(3))
446 cross(ixi^s,3)=two*tmp(ixi^s)
447 v(ixi^s)=w(ixi^s,mom(1))
448 cross(ixi^s,3)=cross(ixi^s,3)+two*v(ixi^s)/x(ixi^s,1)
449 v(ixi^s)=w(ixi^s,mom(2))
450 cross(ixi^s,3)=cross(ixi^s,3)+two*v(ixi^s)/(x(ixi^s,1)*dtan(x(ixi^s,2)))
454 call mpistop(
"Unknown geometry specified")
464 integer,
intent(in) :: ixI^L, ixO^L, idim
465 double precision,
intent(in) :: w(ixI^S, 1:nw), x(ixI^S, 1:^ND)
466 double precision,
intent(out) :: cross(ixI^S,ndir)
468 double precision :: tmp(ixI^S), v(ixI^S)
470 v(ixi^s)=w(ixi^s,mom(idim))
472 call gradient(v,ixi^l,ixo^l,idir,tmp)
473 cross(ixo^s,idir)=tmp(ixo^s)
476 v(ixi^s)=w(ixi^s,mom(idir))
477 call gradient(v,ixi^l,ixo^l,idim,tmp)
478 cross(ixo^s,idir)=cross(ixo^s,idir)+tmp(ixo^s)
487 integer,
intent(in) :: ixI^L, ixO^L
488 double precision,
intent(in) :: qdt, x(ixI^S, 1:ndim)
489 double precision,
intent(inout) :: wCT(ixI^S, 1:nw), w(ixI^S, 1:nw)
493 double precision :: v(ixI^S,1:ndir), vv(ixI^S), divergence(ixI^S)
494 double precision :: tmp(ixI^S),tmp1(ixI^S)
503 v(ixi^s,i)=wct(ixi^s,mom(i))/wct(ixi^s,rho_)
509 call gradient(vv,ixi^l,ixo^l,2,tmp1)
510 tmp(ixo^s)=two*(tmp1(ixo^s)+v(ixi^s,1)/x(ixo^s,1))
513 tmp(ixo^s) = tmp(ixo^s) - (2.d0/3.d0) * divergence(ixo^s)
515 w(ixo^s,mom(1))=w(ixo^s,mom(1))-qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
518 call gradient(vv,ixi^l,ixo^l,2,tmp1)
519 tmp(ixo^s)=tmp1(ixo^s)
520 vv(ixi^s)=v(ixi^s,2)/x(ixi^s,1)
521 call gradient(vv,ixi^l,ixo^l,1,tmp1)
522 tmp(ixo^s)=tmp(ixo^s)+tmp1(ixo^s)*x(ixo^s,1)
524 w(ixo^s,mom(2))=w(ixo^s,mom(2))+qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
529 v(ixi^s,i)=wct(ixi^s,mom(i))/wct(ixi^s,rho_)
535 call gradient(vv,ixi^l,ixo^l,2,tmp1)
536 tmp(ixo^s)=two*(tmp1(ixo^s)+v(ixo^s,1)/x(ixo^s,1))
539 tmp(ixo^s) = tmp(ixo^s) - (2.d0/3.d0) * divergence(ixo^s)
541 w(ixo^s,mom(1))=w(ixo^s,mom(1))-qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
547 call gradient(vv,ixi^l,ixo^l,3,tmp1)
548 tmp(ixo^s)=two*(tmp1(ixo^s)+v(ixo^s,1)/x(ixo^s,1)+v(ixo^s,2)/(x(ixo^s,1)*dtan(x(ixo^s,2))))
551 tmp(ixo^s) = tmp(ixo^s) - (2.d0/3.d0) * divergence(ixo^s)
553 w(ixo^s,mom(1))=w(ixo^s,mom(1))-qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
556 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)))
561 vv(ixi^s)=v(ixi^s,2)/x(ixi^s,1)
562 call gradient(vv,ixi^l,ixo^l,1,tmp1)
563 tmp(ixo^s)=tmp(ixo^s)+tmp1(ixo^s)*x(ixo^s,1)
565 w(ixo^s,mom(2))=w(ixo^s,mom(2))+qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
571 vv(ixi^s)=v(ixi^s,3)/x(ixi^s,1)
572 call gradient(vv,ixi^l,ixo^l,1,tmp1)
573 tmp(ixo^s)=tmp(ixo^s)+tmp1(ixo^s)*x(ixo^s,1)
575 w(ixo^s,mom(3))=w(ixo^s,mom(3))+qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
582 vv(ixi^s)=v(ixi^s,3)/dsin(x(ixi^s,2))
583 call gradient(vv,ixi^l,ixo^l,2,tmp1)
584 tmp(ixo^s)=tmp(ixo^s)+tmp1(ixo^s)*dsin(x(ixo^s,2))
586 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.
integer, 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)