23 call mpistop(
"Hyperdiffusivity only implemented for Cartesian uniform grid")
32 subroutine hyp_coeff(ixI^L, ixO^L, var, idimm, nu_hyp)
34 integer,
intent(in) :: ixI^L, idimm
35 integer,
intent(out) :: ixO^L
36 double precision,
intent(in) :: var(ixI^S)
37 double precision,
intent(out) :: nu_hyp(ixI^S)
39 double precision :: tmp(ixI^S), tmp2(ixI^S)
40 integer :: hx^L, hx1f^L, hx1b^L, hx2b^L
42 double precision,
parameter :: eps=1e-8
45 hxmin^
d=iximin^
d+2*
kr(idimm,^
d);
46 hxmax^
d=iximax^
d-
kr(idimm,^
d);
48 hx1f^l=hx^l+
kr(idimm,^
d);
49 hx1b^l=hx^l-
kr(idimm,^
d);
50 hx2b^l=hx^l-2*
kr(idimm,^
d);
53 tmp(hx^s)= abs(3d0 * (var(hx^s) - var(hx1b^s)) - (var(hx1f^s)-var(hx2b^s)))
55 tmp2(hx^s)= abs((var(hx^s) - var(hx1b^s)))
57 ixomin^
d=hxmin^
d+
kr(idimm,^
d);
58 ixomax^
d=hxmax^
d-
kr(idimm,^
d);
59 hx1f^l=ixo^l+
kr(idimm,^
d);
60 hx1b^l=ixo^l-
kr(idimm,^
d);
62 nu_hyp(ixo^s) = max(tmp(hx1b^s),tmp(ixo^s),tmp(hx1f^s))/max(tmp2(hx1b^s),tmp2(ixo^s),tmp2(hx1f^s),eps)
72 integer,
intent(in) :: ixI^L, idimm
73 integer,
intent(out) :: ixO^L
74 double precision,
intent(in) :: vel(ixI^S,1:ndir)
75 double precision,
intent(out) :: nu_vel(ixI^S)
76 integer :: hx1f^L, hx1b^L
83 ixomin^
d=ixomin^
d+
kr(idimm,^
d);
89 hx1b^l=ixo^l-
kr(ii,^
d);
91 nu_vel(ixo^s) = nu_vel(ixo^s) + (vel(ixo^s,ii) - vel(hx1b^s,ii))/
dxlevel(ii)
93 hx1f^l=ixo^l+
kr(ii,^
d);
94 nu_vel(ixo^s) = nu_vel(ixo^s) + (vel(hx1f^s,ii) - vel(hx1b^s,ii))/(2d0*
dxlevel(ii))
97 where(nu_vel(ixo^s) < 0d0)
98 nu_vel(ixo^s) = -nu_vel(ixo^s)
115 integer,
intent(in) :: ixI^L, idimm
116 integer,
intent(out) :: ixO^L
117 double precision,
intent(in) :: nu_hyper(ixI^S),var(ixI^S)
118 double precision,
intent(out) :: res(ixI^S)
120 integer :: hxf^L, hxb^L
125 hxf^l=ixo^l+
kr(idimm,^
d);
126 hxb^l=ixo^l-
kr(idimm,^
d);
128 res(ixo^s) = 1d0/(
dxlevel(idimm)**2)*&
129 (nu_hyper(hxf^s) * (var(hxf^s)-var(ixo^s))-&
130 nu_hyper(ixo^s) * (var(ixo^s)-var(hxb^s)))
140 integer,
intent(in) :: ixI^L, idimm
141 integer,
intent(out) :: ixO^L
142 double precision,
intent(in) :: nu_hyper(ixI^S),var(ixI^S), var2(ixI^S)
143 double precision,
intent(out) :: res(ixI^S)
145 integer :: hxf^L, hxb^L
150 hxf^l=ixo^l+
kr(idimm,^
d);
151 hxb^l=ixo^l-
kr(idimm,^
d);
153 res(ixo^s) = 1d0/(2d0*
dxlevel(idimm)**2)*&
154 (nu_hyper(hxf^s) * (var2(hxf^s)+var2(ixo^s)) * (var(hxf^s)-var(ixo^s))-&
155 nu_hyper(ixo^s) * (var2(hxb^s)+var2(ixo^s)) * (var(ixo^s)-var(hxb^s)))
166 integer,
intent(in) :: ixI^L, idimm,idimm2
167 integer,
intent(out) :: ixO^L
168 double precision,
intent(in) :: nu_hyper(ixI^S),var(ixI^S)
169 double precision,
intent(out) :: res(ixI^S)
171 integer :: hxfi^L, hxbi^L, hxfifj^L, hxbibj^L, hxfibj^L, hxbifj^L
176 hxfi^l=ixo^l+
kr(idimm2,^
d);
177 hxbi^l=ixo^l-
kr(idimm2,^
d);
179 hxfifj^l=hxfi^l+
kr(idimm,^
d);
180 hxfibj^l=hxfi^l-
kr(idimm,^
d);
182 hxbifj^l=hxbi^l+
kr(idimm,^
d);
183 hxbibj^l=hxbi^l-
kr(idimm,^
d);
186 ((nu_hyper(hxfifj^s) + nu_hyper(hxfi^s)) * (var(hxfifj^s)-var(hxfibj^s))-&
187 (nu_hyper(hxbifj^s) + nu_hyper(hxbi^s)) * (var(hxbifj^s)-var(hxbibj^s)))
201 integer,
intent(in) :: ixI^L, idimm,idimm2
202 integer,
intent(out) :: ixO^L
203 double precision,
intent(in) :: nu_hyper(ixI^S),var(ixI^S),var2(ixI^S)
204 double precision,
intent(out) :: res(ixI^S)
206 integer :: hxfi^L, hxbi^L, hxfifj^L, hxbibj^L, hxfibj^L, hxbifj^L
211 hxfi^l=ixo^l+
kr(idimm2,^
d);
212 hxbi^l=ixo^l-
kr(idimm2,^
d);
214 hxfifj^l=hxfi^l+
kr(idimm,^
d);
215 hxfibj^l=hxfi^l-
kr(idimm,^
d);
217 hxbifj^l=hxbi^l+
kr(idimm,^
d);
218 hxbibj^l=hxbi^l-
kr(idimm,^
d);
221 ((nu_hyper(hxfifj^s) + nu_hyper(hxfi^s)) * var2(hxfi^s) * (var(hxfifj^s)-var(hxfibj^s))-&
222 (nu_hyper(hxbifj^s) + nu_hyper(hxbi^s)) * var2(hxbi^s) * (var(hxbifj^s)-var(hxbibj^s)))
Module with geometry-related routines (e.g., divergence, curl)
integer, parameter cartesian
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer, parameter ndim
Number of spatial dimensions for grid variables.
double precision, dimension(:), allocatable, parameter d
integer nghostcells
Number of ghost cells surrounding a grid.
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
Subroutines for Roe-type Riemann solver for HD.
subroutine second_same_deriv2(ixIL, ixOL, nu_hyper, var2, var, idimm, res)
subroutine second_cross_deriv(ixIL, ixOL, nu_hyper, var, idimm, idimm2, res)
subroutine div_vel_coeff(ixIL, ixOL, vel, idimm, nu_vel)
subroutine hyp_coeff(ixIL, ixOL, var, idimm, nu_hyp)
subroutine second_cross_deriv2(ixIL, ixOL, nu_hyper, var2, var, idimm, idimm2, res)
subroutine second_same_deriv(ixIL, ixOL, nu_hyper, var, idimm, res)
subroutine hyperdiffusivity_init()
This module defines the procedures of a physics module. It contains function pointers for the various...
logical phys_req_diagonal
Whether the physics routines require diagonal ghost cells, for example for computing a curl.