3#include "AMReX_MultiFab.H"
7#include "Numeric/Stencil.H"
17#include "AMReX_TimeIntegrator.H"
19#if AMREX_SPACEDIM == 2
32 BL_PROFILE(
"Integrator::Hydro::Hydro()");
40 pp.
forbid(
"scheme",
"use integration.type instead");
43 pp.
query_default(
"eta_refinement_criterion", value.eta_refinement_criterion , 0.01);
45 pp.
query_default(
"omega_refinement_criterion", value.omega_refinement_criterion, 0.01);
47 pp.
query_default(
"gradu_refinement_criterion", value.gradu_refinement_criterion, 0.01);
49 pp.
query_default(
"p_refinement_criterion", value.p_refinement_criterion, 1e100);
51 pp.
query_default(
"rho_refinement_criterion", value.rho_refinement_criterion, 1e100);
67 pp_forbid(
"velocity.bc",
"--> momentum.bc");
82 pp_forbid(
"roefix",
"--> solver.roe.entropy_fix");
89 value.RegisterNewFab(value.eta_mf, value.eta_bc, 1, nghost,
"eta",
true,
true);
90 value.RegisterNewFab(value.eta_old_mf, value.eta_bc, 1, nghost,
"eta_old",
true,
true);
91 value.RegisterNewFab(value.etadot_mf, value.eta_bc, 1, nghost,
"etadot",
true,
false);
93 value.RegisterNewFab(value.density_mf, value.density_bc, 1, nghost,
"density",
true ,
true);
94 value.RegisterNewFab(value.density_old_mf, value.density_bc, 1, nghost,
"density_old",
false,
true);
96 value.RegisterNewFab(value.energy_mf, value.energy_bc, 1, nghost,
"energy",
true ,
true);
97 value.RegisterNewFab(value.energy_old_mf, value.energy_bc, 1, nghost,
"energy_old" ,
false,
true);
99 value.RegisterNewFab(value.momentum_mf, value.momentum_bc, 2, nghost,
"momentum",
true ,
true, {
"x",
"y"});
100 value.RegisterNewFab(value.momentum_old_mf, value.momentum_bc, 2, nghost,
"momentum_old",
false,
true);
102 value.RegisterNewFab(value.pressure_mf, &value.bc_nothing, 1, nghost,
"pressure",
true,
false);
103 value.RegisterNewFab(value.velocity_mf, &value.bc_nothing, 2, nghost,
"velocity",
true,
false,{
"x",
"y"});
104 value.RegisterNewFab(value.vorticity_mf, &value.bc_nothing, 1, nghost,
"vorticity",
true,
false);
106 value.RegisterNewFab(value.m0_mf, &value.bc_nothing, 1, 0,
"m0",
true,
false);
107 value.RegisterNewFab(value.u0_mf, &value.bc_nothing, 2, 0,
"u0",
true,
false, {
"x",
"y"});
108 value.RegisterNewFab(value.q_mf, &value.bc_nothing, 2, 0,
"q",
true,
false, {
"x",
"y"});
110 value.RegisterNewFab(value.solid.momentum_mf, &value.neumann_bc_D, 2, nghost,
"solid.momentum",
true,
false, {
"x",
"y"});
111 value.RegisterNewFab(value.solid.density_mf, &value.neumann_bc_1, 1, nghost,
"solid.density",
true,
false);
112 value.RegisterNewFab(value.solid.energy_mf, &value.neumann_bc_1, 1, nghost,
"solid.energy",
true,
false);
114 value.RegisterNewFab(value.Source_mf, &value.bc_nothing, 4, 0,
"Source",
true,
false);
117 pp_forbid(
"Velocity.ic.type",
"--> velocity.ic.type");
118 pp_forbid(
"Pressure.ic",
"--> pressure.ic");
119 pp_forbid(
"SolidMomentum.ic",
"--> solid.momentum.ic");
120 pp_forbid(
"SolidDensity.ic.type",
"--> solid.density.ic.type");
121 pp_forbid(
"SolidEnergy.ic.type",
"--> solid.energy.ic.type");
122 pp_forbid(
"Density.ic.type",
"--> density.ic.type");
123 pp_forbid(
"rho_injected.ic.type",
"no longer using rho_injected use m0 instead");
124 pp.
forbid(
"mdot.ic.type",
"replace mdot with u0");
181void Hydro::Initialize(
int lev)
183 BL_PROFILE(
"Integrator::Hydro::Initialize");
185 eta_ic ->Initialize(lev, eta_mf, 0.0);
186 eta_ic ->Initialize(lev, eta_old_mf, 0.0);
187 etadot_mf[lev] ->setVal(0.0);
191 velocity_ic ->Initialize(lev, velocity_mf, 0.0);
192 pressure_ic ->Initialize(lev, pressure_mf, 0.0);
193 density_ic ->Initialize(lev, density_mf, 0.0);
195 density_ic ->Initialize(lev, density_old_mf, 0.0);
198 solid.density_ic ->Initialize(lev, solid.density_mf, 0.0);
199 solid.momentum_ic->Initialize(lev, solid.momentum_mf, 0.0);
200 solid.energy_ic ->Initialize(lev, solid.energy_mf, 0.0);
202 ic_m0 ->Initialize(lev, m0_mf, 0.0);
203 ic_u0 ->Initialize(lev, u0_mf, 0.0);
204 ic_q ->Initialize(lev, q_mf, 0.0);
206 Source_mf[lev] ->setVal(0.0);
211void Hydro::Mix(
int lev)
213 for (amrex::MFIter mfi(*eta_mf[lev],
true); mfi.isValid(); ++mfi)
215 const amrex::Box& bx = mfi.growntilebox();
232 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
234 rho(i, j, k) = eta(i, j, k) * rho(i, j, k) + (1.0 - eta(i, j, k)) * rho_solid(i, j, k);
235 rho_old(i, j, k) = rho(i, j, k);
237 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));
238 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));
239 M_old(i, j, k, 0) = M(i, j, k, 0);
240 M_old(i, j, k, 1) = M(i, j, k, 1);
243 (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)
245 E_solid(i, j, k) * (1.0 - eta(i, j, k));
246 E_old(i, j, k) = E(i, j, k);
256 eta_ic->Initialize(lev, eta_mf, time);
266 if (dynamictimestep.on)
267 Integrator::DynamicTimestep_Update();
272 amrex::ParallelDescriptor::ReduceRealMax(c_max);
273 amrex::ParallelDescriptor::ReduceRealMax(vx_max);
274 amrex::ParallelDescriptor::ReduceRealMax(vy_max);
276 Set::Scalar new_timestep = cfl / ((c_max + vx_max) / DX[0] + (c_max + vy_max) / DX[1]);
280 SetTimestep(new_timestep);
286 std::swap(eta_old_mf, eta_mf);
287 std::swap(density_old_mf[lev], density_mf[lev]);
288 std::swap(momentum_old_mf[lev], momentum_mf[lev]);
289 std::swap(energy_old_mf[lev], energy_mf[lev]);
295 UpdateEta(lev, time);
296 for (amrex::MFIter mfi(*eta_mf[lev],
true); mfi.isValid(); ++mfi)
298 const amrex::Box& bx = mfi.growntilebox();
299 amrex::Array4<const Set::Scalar>
const& eta_new = (*eta_mf[lev]).array(mfi);
300 amrex::Array4<const Set::Scalar>
const& eta = (*eta_old_mf[lev]).array(mfi);
301 amrex::Array4<Set::Scalar>
const& etadot = (*etadot_mf[lev]).array(mfi);
302 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
304 etadot(i, j, k) = (eta_new(i, j, k) - eta(i, j, k)) / dt;
314 amrex::Vector<amrex::MultiFab> solution_new;
315 solution_new.emplace_back(*density_mf[lev].get(),amrex::MakeType::make_alias,0,1);
316 solution_new.emplace_back(*momentum_mf[lev].get(),amrex::MakeType::make_alias,0,2);
317 solution_new.emplace_back(*energy_mf[lev].get(),amrex::MakeType::make_alias,0,1);
320 amrex::Vector<amrex::MultiFab> solution_old;
321 solution_old.emplace_back(*density_old_mf[lev].get(),amrex::MakeType::make_alias,0,1);
322 solution_old.emplace_back(*momentum_old_mf[lev].get(),amrex::MakeType::make_alias,0,2);
323 solution_old.emplace_back(*energy_old_mf[lev].get(),amrex::MakeType::make_alias,0,1);
326 amrex::TimeIntegrator timeintegrator(solution_new, time);
329 timeintegrator.set_rhs([&](amrex::Vector<amrex::MultiFab> & rhs_mf, amrex::Vector<amrex::MultiFab> & solution_mf,
const Set::Scalar time)
332 rhs_mf[0], rhs_mf[1], rhs_mf[2],
333 solution_mf[0],solution_mf[1],solution_mf[2]);
337 timeintegrator.set_post_stage_action([&](amrex::Vector<amrex::MultiFab> & stage_mf,
Set::Scalar time)
339 density_bc->FillBoundary(stage_mf[0],0,1,time,0);
340 stage_mf[0].FillBoundary(
true);
341 momentum_bc->FillBoundary(stage_mf[1],0,2,time,0);
342 stage_mf[1].FillBoundary(
true);
343 energy_bc->FillBoundary(stage_mf[2],0,1,time,0);
344 stage_mf[2].FillBoundary(
true);
348 timeintegrator.advance(solution_old, solution_new, time, dt);
355 Set::Scalar dt_max = std::numeric_limits<Set::Scalar>::max();
356 for (amrex::MFIter mfi(*eta_mf[lev],
false); mfi.isValid(); ++mfi)
358 const amrex::Box& bx = mfi.validbox();
377 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
379 if (eta(i,j,k) < cutoff)
381 rho_new(i,j,k,0) = rho_solid(i,j,k,0);
382 M_new(i,j,k,0) = M_solid(i,j,k,0);
383 M_new(i,j,k,1) = M_solid(i,j,k,1);
384 E_new(i,j,k,0) = E_solid(i,j,k,0);
388 omega(i, j, k) = eta(i, j, k) * (
gradu(1,0) -
gradu(0,1));
390 if (dynamictimestep.on)
392 *dt_max_handle = std::fabs(cfl * DX[0] / (u(i,j,k,0)*eta(i,j,k) + small));
393 *dt_max_handle = std::min(*dt_max_handle, std::fabs(cfl * DX[1] / (u(i,j,k,1)*eta(i,j,k) + small)));
394 *dt_max_handle = std::min(*dt_max_handle, std::fabs(cfl_v * DX[0]*DX[0] / (Source(i,j,k,1)+small)));
395 *dt_max_handle = std::min(*dt_max_handle, std::fabs(cfl_v * DX[1]*DX[1] / (Source(i,j,k,2)+small)));
401 if (dynamictimestep.on)
403 this->DynamicTimestep_SyncTimeStep(lev,dt_max);
410 amrex::MultiFab &rho_rhs_mf,
411 amrex::MultiFab &M_rhs_mf,
412 amrex::MultiFab &E_rhs_mf,
413 const amrex::MultiFab &rho_mf,
414 const amrex::MultiFab &M_mf,
415 const amrex::MultiFab &E_mf)
418 for (amrex::MFIter mfi(*eta_mf[lev],
true); mfi.isValid(); ++mfi)
420 const amrex::Box& bx = mfi.growntilebox();
421 amrex::Array4<const Set::Scalar>
const& eta = (*eta_old_mf[lev]).array(mfi);
434 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
437 Set::Scalar etarho_fluid = rho(i,j,k) - (1.-eta(i,j,k)) * rho_solid(i,j,k);
438 Set::Scalar etaE_fluid = E(i,j,k) - (1.-eta(i,j,k)) * E_solid(i,j,k);
440 Set::Vector etaM_fluid( M(i,j,k,0) - (1.-eta(i,j,k)) * M_solid(i,j,k,0),
441 M(i,j,k,1) - (1.-eta(i,j,k)) * M_solid(i,j,k,1) );
445 v(i,j,k,0) = etaM_fluid(0) / (etarho_fluid + small);
446 v(i,j,k,1) = etaM_fluid(1) / (etarho_fluid + small);
447 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;
452 amrex::Box domain = geom[lev].Domain();
454 for (amrex::MFIter mfi(*eta_mf[lev],
false); mfi.isValid(); ++mfi)
456 const amrex::Box& bx = mfi.validbox();
487 amrex::Array4<Set::Scalar>
const& Source = (*Source_mf[lev]).array(mfi);
489 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
516 for (
int p = 0; p < 2; p++)
517 for (
int q = 0; q < 2; q++)
518 for (
int r = 0; r < 2; r++)
521 (hess_M(r,p,q) -
gradu(r,q)*gradrho(p) -
gradu(r,p)*gradrho(q) - u(r)*hess_rho(p,q))
527 for (
int p = 0; p<2; p++)
528 for (
int q = 0; q<2; q++)
529 for (
int r = 0; r<2; r++)
530 for (
int s = 0; s<2; s++)
533 if (p==r && q==s) Mpqrs += 0.5 * mu;
535 Ldot0(p) += 0.5*Mpqrs * (u(r) - u0(r)) * hess_eta(q, s);
536 div_tau(p) += 2.0*Mpqrs * hess_u(r,s,q);
539 Source(i,j, k, 0) = mdot0;
540 Source(i,j, k, 1) = Pdot0(0) - Ldot0(0);
541 Source(i,j, k, 2) = Pdot0(1) - Ldot0(1);
542 Source(i,j, k, 3) = qdot0;
545 Source(i,j,k,1) -= lagrange*(u-u0).dot(grad_eta)*grad_eta(0);
546 Source(i,j,k,2) -= lagrange*(u-u0).dot(grad_eta)*grad_eta(1);
550 const int X = 0, Y = 1;
581 flux_xlo = riemannsolver->Solve(state_xlo_fluid, state_x_fluid, gamma, pref, small) * eta(i,j,k);
582 flux_ylo = riemannsolver->Solve(state_ylo_fluid, state_y_fluid, gamma, pref, small) * eta(i,j,k);
585 flux_xhi = riemannsolver->Solve(state_x_fluid, state_xhi_fluid, gamma, pref, small) * eta(i,j,k);
586 flux_yhi = riemannsolver->Solve(state_y_fluid, state_yhi_fluid, gamma, pref, small) * eta(i,j,k);
597 (flux_xlo.
mass - flux_xhi.
mass) / DX[0] +
598 (flux_ylo.
mass - flux_yhi.
mass) / DX[1] +
606 etadot(i,j,k) * (rho(i,j,k) - rho_solid(i,j,k)) / (eta(i,j,k) + small)
615 div_tau(0) * eta(i,j,k) +
623 etadot(i,j,k)*(M(i,j,k,0) - M_solid(i,j,k,0)) / (eta(i,j,k) + small)
630 div_tau(1) * eta(i,j,k) +
638 etadot(i,j,k)*(M(i,j,k,1) - M_solid(i,j,k,1)) / (eta(i,j,k)+small)
652 etadot(i,j,k)*(E(i,j,k) - E_solid(i,j,k)) / (eta(i,j,k)+small)
657 if ((rho_rhs(i,j,k) != rho_rhs(i,j,k)) ||
658 (M_rhs(i,j,k,0) != M_rhs(i,j,k,0)) ||
659 (M_rhs(i,j,k,1) != M_rhs(i,j,k,1)) ||
660 (E_rhs(i,j,k) != E_rhs(i,j,k)))
719 omega(i, j, k) = eta(i, j, k) * (
gradu(1,0) -
gradu(0,1));
726 BL_PROFILE(
"Integrator::Hydro::Regrid");
727 Source_mf[lev]->setVal(0.0);
728 if (lev < finest_level)
return;
734void Hydro::TagCellsForRefinement(
int lev, amrex::TagBoxArray& a_tags,
Set::Scalar,
int)
736 BL_PROFILE(
"Integrator::Flame::TagCellsForRefinement");
739 Set::Scalar dr = sqrt(AMREX_D_TERM(DX[0] * DX[0], +DX[1] * DX[1], +DX[2] * DX[2]));
742 for (amrex::MFIter mfi(*eta_mf[lev],
true); mfi.isValid(); ++mfi) {
743 const amrex::Box& bx = mfi.tilebox();
744 amrex::Array4<char>
const& tags = a_tags.array(mfi);
745 amrex::Array4<const Set::Scalar>
const& eta = (*eta_mf[lev]).array(mfi);
747 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
749 if (grad_eta.lpNorm<2>() * dr * 2 > eta_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
754 for (amrex::MFIter mfi(*vorticity_mf[lev],
true); mfi.isValid(); ++mfi) {
755 const amrex::Box& bx = mfi.tilebox();
756 amrex::Array4<char>
const& tags = a_tags.array(mfi);
757 amrex::Array4<const Set::Scalar>
const& omega = (*vorticity_mf[lev]).array(mfi);
759 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
762 if (grad_omega.lpNorm<2>() * dr * 2 > omega_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
767 for (amrex::MFIter mfi(*velocity_mf[lev],
true); mfi.isValid(); ++mfi) {
768 const amrex::Box& bx = mfi.tilebox();
769 amrex::Array4<char>
const& tags = a_tags.array(mfi);
770 amrex::Array4<const Set::Scalar>
const& v = (*velocity_mf[lev]).array(mfi);
772 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
775 if (grad_u.lpNorm<2>() * dr * 2 > gradu_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
780 for (amrex::MFIter mfi(*pressure_mf[lev],
true); mfi.isValid(); ++mfi) {
781 const amrex::Box& bx = mfi.tilebox();
782 amrex::Array4<char>
const& tags = a_tags.array(mfi);
783 amrex::Array4<const Set::Scalar>
const& p = (*pressure_mf[lev]).array(mfi);
785 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
788 if (grad_p.lpNorm<2>() * dr * 2 > p_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
793 for (amrex::MFIter mfi(*density_mf[lev],
true); mfi.isValid(); ++mfi) {
794 const amrex::Box& bx = mfi.tilebox();
795 amrex::Array4<char>
const& tags = a_tags.array(mfi);
796 amrex::Array4<const Set::Scalar>
const& rho = (*density_mf[lev]).array(mfi);
798 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
801 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(...)
#define pp_queryclass(...)
Initialize Laminates in a matrix.
void forbid(std::string name, std::string explanation)
int AnyUnusedInputs(bool inscopeonly=true, bool verbose=false)
void select_default(std::string name, PTRTYPE *&ic_eta, Args &&... args)
static int AllUnusedInputs()
int query_default(std::string name, T &value, T defaultvalue)
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 Warning(std::string file, std::string func, int line, 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