6#include "Numeric/Stencil.H"
14#if AMREX_SPACEDIM == 2
27 BL_PROFILE(
"Integrator::Hydro::Hydro()");
36 pp.
query_default(
"eta_refinement_criterion", value.eta_refinement_criterion , 0.01);
38 pp.
query_default(
"omega_refinement_criterion", value.omega_refinement_criterion, 0.01);
40 pp.
query_default(
"gradu_refinement_criterion", value.gradu_refinement_criterion, 0.01);
42 pp.
query_default(
"p_refinement_criterion", value.p_refinement_criterion, 1e100);
44 pp.
query_default(
"rho_refinement_criterion", value.rho_refinement_criterion, 1e100);
60 pp_forbid(
"velocity.bc",
"--> momentum.bc");
75 pp_forbid(
"roefix",
"--> solver.roe.entropy_fix");
82 value.RegisterNewFab(value.eta_mf, value.eta_bc, 1, nghost,
"eta",
true );
83 value.RegisterNewFab(value.eta_old_mf, value.eta_bc, 1, nghost,
"eta_old",
true);
84 value.RegisterNewFab(value.etadot_mf, value.eta_bc, 1, nghost,
"etadot",
true );
86 value.RegisterNewFab(value.density_mf, value.density_bc, 1, nghost,
"density",
true );
87 value.RegisterNewFab(value.density_old_mf, value.density_bc, 1, nghost,
"density_old",
false);
89 value.RegisterNewFab(value.energy_mf, value.energy_bc, 1, nghost,
"energy",
true );
90 value.RegisterNewFab(value.energy_old_mf, value.energy_bc, 1, nghost,
"energy_old" ,
false);
92 value.RegisterNewFab(value.momentum_mf, value.momentum_bc, 2, nghost,
"momentum",
true , {
"x",
"y"});
93 value.RegisterNewFab(value.momentum_old_mf, value.momentum_bc, 2, nghost,
"momentum_old",
false);
95 value.RegisterNewFab(value.pressure_mf, &value.bc_nothing, 1, nghost,
"pressure",
true);
96 value.RegisterNewFab(value.velocity_mf, &value.bc_nothing, 2, nghost,
"velocity",
true, {
"x",
"y"});
97 value.RegisterNewFab(value.vorticity_mf, &value.bc_nothing, 1, nghost,
"vorticity",
true);
99 value.RegisterNewFab(value.m0_mf, &value.bc_nothing, 1, 0,
"m0",
true);
100 value.RegisterNewFab(value.u0_mf, &value.bc_nothing, 2, 0,
"u0",
true, {
"x",
"y"});
101 value.RegisterNewFab(value.q_mf, &value.bc_nothing, 2, 0,
"q",
true, {
"x",
"y"});
103 value.RegisterNewFab(value.solid.momentum_mf, &value.neumann_bc_D, 2, nghost,
"solid.momentum",
true, {
"x",
"y"});
104 value.RegisterNewFab(value.solid.density_mf, &value.neumann_bc_1, 1, nghost,
"solid.density",
true);
105 value.RegisterNewFab(value.solid.energy_mf, &value.neumann_bc_1, 1, nghost,
"solid.energy",
true);
107 value.RegisterNewFab(value.Source_mf, &value.bc_nothing, 4, 0,
"Source",
true);
110 pp_forbid(
"Velocity.ic.type",
"--> velocity.ic.type");
111 pp_forbid(
"Pressure.ic",
"--> pressure.ic");
112 pp_forbid(
"SolidMomentum.ic",
"--> solid.momentum.ic");
113 pp_forbid(
"SolidDensity.ic.type",
"--> solid.density.ic.type");
114 pp_forbid(
"SolidEnergy.ic.type",
"--> solid.energy.ic.type");
115 pp_forbid(
"Density.ic.type",
"--> density.ic.type");
116 pp_forbid(
"rho_injected.ic.type",
"no longer using rho_injected use m0 instead");
117 pp.
forbid(
"mdot.ic.type",
"replace mdot with u0");
159void Hydro::Initialize(
int lev)
161 BL_PROFILE(
"Integrator::Hydro::Initialize");
163 eta_ic ->Initialize(lev, eta_mf, 0.0);
164 eta_ic ->Initialize(lev, eta_old_mf, 0.0);
165 etadot_mf[lev] ->setVal(0.0);
169 velocity_ic ->Initialize(lev, velocity_mf, 0.0);
170 pressure_ic ->Initialize(lev, pressure_mf, 0.0);
171 density_ic ->Initialize(lev, density_mf, 0.0);
173 density_ic ->Initialize(lev, density_old_mf, 0.0);
176 solid.density_ic ->Initialize(lev, solid.density_mf, 0.0);
177 solid.momentum_ic->Initialize(lev, solid.momentum_mf, 0.0);
178 solid.energy_ic ->Initialize(lev, solid.energy_mf, 0.0);
180 ic_m0 ->Initialize(lev, m0_mf, 0.0);
181 ic_u0 ->Initialize(lev, u0_mf, 0.0);
182 ic_q ->Initialize(lev, q_mf, 0.0);
184 Source_mf[lev] ->setVal(0.0);
189void Hydro::Mix(
int lev)
191 for (amrex::MFIter mfi(*eta_mf[lev],
true); mfi.isValid(); ++mfi)
193 const amrex::Box& bx = mfi.growntilebox();
210 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
212 rho(i, j, k) = eta(i, j, k) * rho(i, j, k) + (1.0 - eta(i, j, k)) * rho_solid(i, j, k);
213 rho_old(i, j, k) = rho(i, j, k);
215 M(i, j, k, 0) = (rho(i, j, k)*v(i, j, k, 0))*eta(i, j, k) + M_solid(i, j, k, 0)*(1.0-eta(i, j, k));
216 M(i, j, k, 1) = (rho(i, j, k)*v(i, j, k, 1))*eta(i, j, k) + M_solid(i, j, k, 1)*(1.0-eta(i, j, k));
217 M_old(i, j, k, 0) = M(i, j, k, 0);
218 M_old(i, j, k, 1) = M(i, j, k, 1);
221 (0.5 * (v(i, j, k, 0) * v(i, j, k, 0) + v(i, j, k, 1) * v(i, j, k, 1)) * rho(i, j, k) + p(i, j, k) / (gamma - 1.0)) * eta(i, j, k)
223 E_solid(i, j, k) * (1.0 - eta(i, j, k));
224 E_old(i, j, k) = E(i, j, k);
234 eta_ic->Initialize(lev, eta_mf, time);
244 Integrator::DynamicTimestep_Update();
250 amrex::ParallelDescriptor::ReduceRealMax(c_max);
251 amrex::ParallelDescriptor::ReduceRealMax(vx_max);
252 amrex::ParallelDescriptor::ReduceRealMax(vy_max);
254 Set::Scalar new_timestep = cfl / ((c_max + vx_max) / DX[0] + (c_max + vy_max) / DX[1]);
258 SetTimestep(new_timestep);
264 std::swap(eta_old_mf, eta_mf);
265 std::swap(density_old_mf[lev], density_mf[lev]);
266 std::swap(momentum_old_mf[lev], momentum_mf[lev]);
267 std::swap(energy_old_mf[lev], energy_mf[lev]);
268 Set::Scalar dt_max = std::numeric_limits<Set::Scalar>::max();
270 UpdateEta(lev, time);
272 for (amrex::MFIter mfi(*eta_mf[lev],
true); mfi.isValid(); ++mfi)
274 const amrex::Box& bx = mfi.growntilebox();
275 amrex::Array4<const Set::Scalar>
const& eta_new = (*eta_mf[lev]).array(mfi);
276 amrex::Array4<const Set::Scalar>
const& eta = (*eta_old_mf[lev]).array(mfi);
277 amrex::Array4<Set::Scalar>
const& etadot = (*etadot_mf[lev]).array(mfi);
291 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
293 etadot(i, j, k) = (eta_new(i, j, k) - eta(i, j, k)) / dt;
295 Set::Scalar etarho_fluid = rho(i,j,k) - (1.-eta(i,j,k)) * rho_solid(i,j,k);
296 Set::Scalar etaE_fluid = E(i,j,k) - (1.-eta(i,j,k)) * E_solid(i,j,k);
298 Set::Vector etaM_fluid( M(i,j,k,0) - (1.-eta(i,j,k)) * M_solid(i,j,k,0),
299 M(i,j,k,1) - (1.-eta(i,j,k)) * M_solid(i,j,k,1) );
303 v(i,j,k,0) = etaM_fluid(0) / (etarho_fluid + small);
304 v(i,j,k,1) = etaM_fluid(1) / (etarho_fluid + small);
306 p(i,j,k) = (etaE_fluid / (eta(i, j, k) + small) - 0.5 * (etaM_fluid(0)*etaM_fluid(0) + etaM_fluid(1)*etaM_fluid(1)) / (etarho_fluid + small)) * ((gamma - 1.0) / (eta(i, j, k) + small))-pref;
311 amrex::Box domain = geom[lev].Domain();
313 for (amrex::MFIter mfi(*eta_mf[lev],
false); mfi.isValid(); ++mfi)
315 const amrex::Box& bx = mfi.validbox();
339 amrex::Array4<Set::Scalar>
const& Source = (*Source_mf[lev]).array(mfi);
343 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
369 for (
int p = 0; p < 2; p++)
370 for (
int q = 0; q < 2; q++)
371 for (
int r = 0; r < 2; r++)
374 (hess_M(r,p,q) -
gradu(r,q)*gradrho(p) -
gradu(r,p)*gradrho(q) - u(r)*hess_rho(p,q))
380 for (
int p = 0; p<2; p++)
381 for (
int q = 0; q<2; q++)
382 for (
int r = 0; r<2; r++)
383 for (
int s = 0; s<2; s++)
386 if (p==r && q==s) Mpqrs += 0.5 * mu;
388 Ldot0(p) += 0.5*Mpqrs * (u(r) - u0(r)) * hess_eta(q, s);
389 div_tau(p) += 2.0*Mpqrs * hess_u(r,s,q);
392 Source(i,j, k, 0) = mdot0;
393 Source(i,j, k, 1) = Pdot0(0) - Ldot0(0);
394 Source(i,j, k, 2) = Pdot0(1) - Ldot0(1);
395 Source(i,j, k, 3) = qdot0;
398 Source(i,j,k,1) -= lagrange*u.dot(grad_eta)*grad_eta(0);
399 Source(i,j,k,2) -= lagrange*u.dot(grad_eta)*grad_eta(1);
403 const int X = 0, Y = 1;
434 flux_xlo = roesolver->Solve(state_xlo_fluid, state_x_fluid, gamma, pref, small) * eta(i,j,k);
435 flux_ylo = roesolver->Solve(state_ylo_fluid, state_y_fluid, gamma, pref, small) * eta(i,j,k);
438 flux_xhi = roesolver->Solve(state_x_fluid, state_xhi_fluid, gamma, pref, small) * eta(i,j,k);
439 flux_yhi = roesolver->Solve(state_y_fluid, state_yhi_fluid, gamma, pref, small) * eta(i,j,k);
450 (flux_xlo.
mass - flux_xhi.
mass) / DX[0] +
451 (flux_ylo.
mass - flux_yhi.
mass) / DX[1] +
454 rho_new(i, j, k) = rho(i, j, k) +
458 etadot(i,j,k) * (rho(i,j,k) - rho_solid(i,j,k)) / (eta(i,j,k) + small)
461 if (rho_new(i,j,k) != rho_new(i,j,k))
491 div_tau(0) * eta(i,j,k) +
495 M_new(i, j, k, 0) = M(i, j, k, 0) +
499 etadot(i,j,k)*(M(i,j,k,0) - M_solid(i,j,k,0)) / (eta(i,j,k) + small)
505 div_tau(1) * eta(i,j,k) +
509 M_new(i, j, k, 1) = M(i, j, k, 1) +
513 etadot(i,j,k)*(M(i,j,k,1) - M_solid(i,j,k,1)) / (eta(i,j,k)+small)
521 E_new(i, j, k) = E(i, j, k) +
525 etadot(i,j,k)*(E(i,j,k) - E_solid(i,j,k)) / (eta(i,j,k)+small)
529 if (eta(i,j,k) < cutoff)
531 rho_new(i,j,k,0) = rho_solid(i,j,k,0);
532 M_new(i,j,k,0) = M_solid(i,j,k,0);
533 M_new(i,j,k,1) = M_solid(i,j,k,1);
534 E_new(i,j,k,0) = E_solid(i,j,k,0);
541 *dt_max_handle = std::fabs(cfl * DX[0] / (u(0)*eta(i,j,k) + small));
542 *dt_max_handle = std::min(*dt_max_handle, std::fabs(cfl * DX[1] / (u(1)*eta(i,j,k) + small)));
543 *dt_max_handle = std::min(*dt_max_handle, std::fabs(cfl_v * DX[0]*DX[0] / (Source(i,j,k,1)+small)));
544 *dt_max_handle = std::min(*dt_max_handle, std::fabs(cfl_v * DX[1]*DX[1] / (Source(i,j,k,2)+small)));
547 omega(i, j, k) = eta(i, j, k) * (
gradu(1,0) -
gradu(0,1));
552 this->DynamicTimestep_SyncTimeStep(lev,dt_max);
558 BL_PROFILE(
"Integrator::Hydro::Regrid");
559 Source_mf[lev]->setVal(0.0);
560 if (lev < finest_level)
return;
566void Hydro::TagCellsForRefinement(
int lev, amrex::TagBoxArray& a_tags,
Set::Scalar,
int)
568 BL_PROFILE(
"Integrator::Flame::TagCellsForRefinement");
571 Set::Scalar dr = sqrt(AMREX_D_TERM(DX[0] * DX[0], +DX[1] * DX[1], +DX[2] * DX[2]));
574 for (amrex::MFIter mfi(*eta_mf[lev],
true); mfi.isValid(); ++mfi) {
575 const amrex::Box& bx = mfi.tilebox();
576 amrex::Array4<char>
const& tags = a_tags.array(mfi);
577 amrex::Array4<const Set::Scalar>
const& eta = (*eta_mf[lev]).array(mfi);
579 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
581 if (grad_eta.lpNorm<2>() * dr * 2 > eta_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
586 for (amrex::MFIter mfi(*vorticity_mf[lev],
true); mfi.isValid(); ++mfi) {
587 const amrex::Box& bx = mfi.tilebox();
588 amrex::Array4<char>
const& tags = a_tags.array(mfi);
589 amrex::Array4<const Set::Scalar>
const& omega = (*vorticity_mf[lev]).array(mfi);
591 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
594 if (grad_omega.lpNorm<2>() * dr * 2 > omega_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
599 for (amrex::MFIter mfi(*velocity_mf[lev],
true); mfi.isValid(); ++mfi) {
600 const amrex::Box& bx = mfi.tilebox();
601 amrex::Array4<char>
const& tags = a_tags.array(mfi);
602 amrex::Array4<const Set::Scalar>
const& v = (*velocity_mf[lev]).array(mfi);
604 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
607 if (grad_u.lpNorm<2>() * dr * 2 > gradu_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
612 for (amrex::MFIter mfi(*pressure_mf[lev],
true); mfi.isValid(); ++mfi) {
613 const amrex::Box& bx = mfi.tilebox();
614 amrex::Array4<char>
const& tags = a_tags.array(mfi);
615 amrex::Array4<const Set::Scalar>
const& p = (*pressure_mf[lev]).array(mfi);
617 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
620 if (grad_p.lpNorm<2>() * dr * 2 > p_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
625 for (amrex::MFIter mfi(*density_mf[lev],
true); mfi.isValid(); ++mfi) {
626 const amrex::Box& bx = mfi.tilebox();
627 amrex::Array4<char>
const& tags = a_tags.array(mfi);
628 amrex::Array4<const Set::Scalar>
const& rho = (*density_mf[lev]).array(mfi);
630 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
633 if (grad_rho.lpNorm<2>() * dr * 2 > rho_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
#define pp_query_required(...)
#define pp_query_default(...)
Initialize Laminates in a matrix.
void forbid(std::string name, std::string explanation, std::string file="", std::string func="", int line=-1)
void queryclass(std::string name, T *value, std::string file="", std::string func="", int line=-1)
void select_default(std::string name, PTRTYPE *&ic_eta, Args &&... args)
int query_default(std::string name, T &value, T defaultvalue, std::string="", std::string="", int=-1)
Roe Riemann Solver based on Gas Dynamics - Culbert B. Laney.
Collection of numerical integrator objects.
static AMREX_FORCE_INLINE std::array< StencilType, AMREX_SPACEDIM > GetStencil(const int i, const int j, const int k, const amrex::Box domain)
AMREX_FORCE_INLINE Set::Matrix Hessian(const amrex::Array4< const Set::Scalar > &f, const int &i, const int &j, const int &k, const int &m, const Set::Scalar dx[AMREX_SPACEDIM], std::array< StencilType, AMREX_SPACEDIM > stencil=DefaultType)
AMREX_FORCE_INLINE Set::Vector Gradient(const amrex::Array4< const Set::Scalar > &f, const int &i, const int &j, const int &k, const int &m, const Set::Scalar dx[AMREX_SPACEDIM], std::array< StencilType, AMREX_SPACEDIM > stencil=DefaultType)
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
void ParallelMessage(std::string file, std::string func, int line, Args const &... args)
void Abort(const char *msg)
AMREX_FORCE_INLINE void Assert(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)
void Message(std::string file, std::string func, int line, Args const &... args)
void Exception(std::string file, std::string func, int line, Args const &... args)
Set::Scalar momentum_normal
Set::Scalar momentum_tangent