Alamo
Hydro.cpp
Go to the documentation of this file.
1
2#include "Hydro.H"
3#include "AMReX_MultiFab.H"
4#include "IO/ParmParse.H"
5#include "BC/Constant.H"
6#include "BC/Expression.H"
7#include "Numeric/Stencil.H"
8#include "IC/Constant.H"
9#include "IC/Laminate.H"
10#include "IC/Expression.H"
11#include "IC/BMP.H"
12#include "IC/PNG.H"
16#include "AMReX_TimeIntegrator.H"
17
18#include "Model/Gas/Gas.H"
23#include "Model/Gas/EOS/EOS.H"
24#include "Model/Gas/EOS/CPG.H"
25
26namespace Integrator
27{
28
30{
31 pp_queryclass(*this);
32}
33
34void
36{
37 BL_PROFILE("Integrator::Hydro::Hydro()");
38 {
39 // pp.query_default("r_refinement_criterion", value.r_refinement_criterion , 0.01);
40 // energy-based refinement
41 // pp.query_default("e_refinement_criterion", value.e_refinement_criterion , 0.01);
42 // momentum-based refinement
43 // pp.query_default("m_refinement_criterion", value.m_refinement_criterion , 0.01);
44
45 pp.forbid("scheme","use integration.type instead");
46
47 // eta-based refinement
48 pp.query_default("eta_refinement_criterion", value.eta_refinement_criterion , 0.01);
49 // vorticity-based refinement
50 pp.query_default("omega_refinement_criterion", value.omega_refinement_criterion, 0.01);
51 // velocity gradient-based refinement
52 pp.query_default("gradu_refinement_criterion", value.gradu_refinement_criterion, 0.01);
53 // pressure-based refinement
54 pp.query_default("p_refinement_criterion", value.p_refinement_criterion, 1e100);
55 // density-based refinement
56 pp.query_default("rho_refinement_criterion", value.rho_refinement_criterion, 1e100);
57
58 pp_forbid("gamma", "replaced by gas->gamma(...)"); // gamma for gamma law
59 pp_query_required("cfl", value.cfl); // cfl condition
60 pp_query_default("cfl_v", value.cfl_v,1E100); // cfl condition
61 pp_forbid("mu", "replaced with gas->dynamic_viscosity(...)"); // linear viscosity coefficient
62 pp_forbid("Lfactor","replaced with mu");
63 //pp_query_default("Lfactor", value.Lfactor,1.0); // (to be removed) test factor for viscous source
64 pp_forbid("Pfactor","replaced with mu");
65 //pp_query_default("Pfactor", value.Pfactor,1.0); // (to be removed) test factor for viscous source
66 pp_forbid("pref", "deprecated - use absolute pressure"); // reference pressure for Roe solver
67
68 pp_forbid("rho.bc","--> density.bc");
69 pp_forbid("p.bc","--> pressure.bc");
70 pp_forbid("v.bc", "--> velocity.bc");
71 pp_forbid("pressure.bc","--> energy.bc");
72 pp_forbid("velocity.bc","--> momentum.bc");
73
74 // Boundary condition for density
75 pp.select_default<BC::Constant,BC::Expression>("density.bc",value.density_bc,1);
76 // Boundary condition for energy
77 pp.select_default<BC::Constant,BC::Expression>("energy.bc",value.energy_bc,1);
78 // Boundary condition for momentum
79 pp.select_default<BC::Constant,BC::Expression>("momentum.bc",value.momentum_bc,2);
80
81 if (!value.managed)
82 {
83 // Boundary condition for phase field order parameter
84 pp.select_default<BC::Constant,BC::Expression>("pf.eta.bc",value.eta_bc,1);
85 }
86
87 pp_query_default("small",value.small,1E-8); // small regularization value
88 pp_query_default("cutoff",value.cutoff,-1E100); // cutoff value
89 pp_query_default("lagrange",value.lagrange,0.0); // lagrange no-penetration factor
90
91 pp_forbid("roefix","--> solver.roe.entropy_fix"); // Roe solver entropy fix
92
93 }
94 // Register FabFields:
95 {
96 int nghost = 1;
97
98 if (!value.managed)
99 {
100 value.eta_mf = new Set::Field<Set::Scalar>();
102 value.RegisterNewFab(*value.eta_mf, value.eta_bc, 1, nghost, "eta", true, true);
103 value.RegisterNewFab(*value.eta_old_mf, value.eta_bc, 1, nghost, "eta_old", true, true);
104 }
105 value.RegisterNewFab(value.etadot_mf, value.eta_bc, 1, nghost, "etadot", true, false);
106
107 value.RegisterNewFab(value.density_mf, value.density_bc, 1, nghost, "density", true , true);
108 value.RegisterNewFab(value.density_old_mf, value.density_bc, 1, nghost, "density_old", false, true);
109
110 value.RegisterNewFab(value.energy_mf, value.energy_bc, 1, nghost, "energy", true ,true);
111 value.RegisterNewFab(value.energy_old_mf, value.energy_bc, 1, nghost, "energy_old" , false, true);
112
113 value.RegisterNewFab(value.momentum_mf, value.momentum_bc, 2, nghost, "momentum", true ,true, {"x","y"});
114 value.RegisterNewFab(value.momentum_old_mf, value.momentum_bc, 2, nghost, "momentum_old", false, true);
115
116 value.RegisterNewFab(value.pressure_mf, &value.bc_nothing, 1, nghost, "pressure", true, false);
117 value.RegisterNewFab(value.temperature_mf, &value.bc_nothing, 1, nghost, "temperature", true, false);
118 value.RegisterNewFab(value.velocity_mf, &value.bc_nothing, 2, nghost, "velocity", true, false,{"x","y"});
119 value.RegisterNewFab(value.vorticity_mf, &value.bc_nothing, 1, nghost, "vorticity", true, false);
120
121 value.RegisterNewFab(value.m0_mf, &value.bc_nothing, 1, 0, "m0", true, false);
122 value.RegisterNewFab(value.u0_mf, &value.bc_nothing, 2, 0, "u0", true, false, {"x","y"});
123 value.RegisterNewFab(value.q_mf, &value.bc_nothing, 2, 0, "q", true, false, {"x","y"});
124
125 value.RegisterNewFab(value.solid.momentum_mf, &value.neumann_bc_D, 2, nghost, "solid.momentum", true, false, {"x","y"});
126 value.RegisterNewFab(value.solid.density_mf, &value.neumann_bc_1, 1, nghost, "solid.density", true, false);
127 value.RegisterNewFab(value.solid.energy_mf, &value.neumann_bc_1, 1, nghost, "solid.energy", true, false);
128
129 value.RegisterNewFab(value.Source_mf, &value.bc_nothing, 4, 0, "Source", true, false);
130
131 value.RegisterNewFab(value.mass_fraction_mf, &value.bc_nothing, 1, nghost, "mass_fraction", true , true);
132 value.RegisterNewFab(value.mole_fraction_mf, &value.bc_nothing, 1, nghost, "mole_fraction", true , true);
133 value.RegisterNewFab(value.scratch_mf, &value.bc_nothing, 1, nghost, "scratch", false , false);
134 }
135
136 pp_forbid("Velocity.ic.type", "--> velocity.ic.type");
137 pp_forbid("Pressure.ic", "--> pressure.ic");
138 pp_forbid("SolidMomentum.ic", "--> solid.momentum.ic");
139 pp_forbid("SolidDensity.ic.type", "--> solid.density.ic.type");
140 pp_forbid("SolidEnergy.ic.type", "--> solid.energy.ic.type");
141 pp_forbid("Density.ic.type", "--> density.ic.type");
142 pp_forbid("rho_injected.ic.type","no longer using rho_injected use m0 instead");
143 pp.forbid("mdot.ic.type", "replace mdot with u0");
144
145
146 // ORDER PARAMETER
147
148 if (!value.managed)
149 {
150 // eta initial condition
152 }
153
154 // PRIMITIVE FIELD INITIAL CONDITIONS
155
156 // velocity initial condition
157 pp.select_default<IC::Constant,IC::Expression>("velocity.ic",value.velocity_ic,value.geom);
158 // solid pressure initial condition
159 pp.select_default<IC::Constant,IC::Expression>("pressure.ic",value.pressure_ic,value.geom);
160 // density initial condition type
161 pp.select_default<IC::Constant,IC::Expression>("density.ic",value.density_ic,value.geom);
162
163
164 // SOLID FIELDS
165
166 // solid momentum initial condition
167 pp.select_default<IC::Constant,IC::Expression>("solid.momentum.ic",value.solid.momentum_ic,value.geom);
168 // solid density initial condition
169 pp.select_default<IC::Constant,IC::Expression>("solid.density.ic",value.solid.density_ic,value.geom);
170 // solid energy initial condition
171 pp.select_default<IC::Constant,IC::Expression>("solid.energy.ic",value.solid.energy_ic,value.geom);
172
173
174 // DIFFUSE BOUNDARY SOURCES
175
176 // diffuse boundary prescribed mass flux
177 pp.select_default<IC::Constant,IC::Expression>("m0.ic",value.ic_m0,value.geom);
178 // diffuse boundary prescribed velocity
179 pp.select_default<IC::Constant,IC::Expression>("u0.ic",value.ic_u0,value.geom);
180 // diffuse boundary prescribed heat flux
181 pp.select_default<IC::Constant,IC::Expression>("q.ic",value.ic_q,value.geom);
182
183 // Riemann solver
187
188 // Gas model (Thermo, Transport, and EOS)
189 pp.queryclass<Model::Gas::Gas>("gas", value.gas);
190 value.nspecies = value.gas.nspecies;
191 std::cout << value.gas.thermo->model_name() << "\n";
192 std::cout << value.gas.transport->model_name() << "\n";
193 std::cout << value.gas.eos->model_name() << "\n";
194 std::cout << value.nspecies << "\n";
195
196 std::string prescribedflowmode_str;
197 //
198 pp.query_validate("prescribedflowmode",prescribedflowmode_str,{"absolute","relative"});
199 if (prescribedflowmode_str == "absolute") value.prescribedflowmode = PrescribedFlowMode::Absolute;
200 else if (prescribedflowmode_str == "relative") value.prescribedflowmode = PrescribedFlowMode::Relative;
201
202 pp.queryarr_default("g",value.g,Set::Vector::Zero());
203
204 bool allow_unused;
205 // Set this to true to allow unused inputs without error.
206 // (Not recommended.)
207 pp.query_default("allow_unused",allow_unused,false);
208 if (!allow_unused && pp.AnyUnusedInputs(true, false))
209 {
210 Util::Warning(INFO,"The following inputs were specified but not used:");
211 pp.AllUnusedInputs();
212 Util::Exception(INFO,"Aborting. Specify 'allow_unused=True` to ignore this error.");
213 }
214}
215
216
217void Hydro::Initialize(int lev)
218{
219 BL_PROFILE("Integrator::Hydro::Initialize");
220
221 if (!managed)
222 {
223 eta_ic ->Initialize(lev, *eta_mf, 0.0);
224 eta_ic ->Initialize(lev, *eta_old_mf, 0.0);
225 }
226 etadot_mf[lev] ->setVal(0.0);
227
228 //flux_mf[lev] ->setVal(0.0);
229
232 density_ic ->Initialize(lev, density_mf, 0.0);
233
235
236 solid.density_ic ->Initialize(lev, solid.density_mf, 0.0);
237 solid.momentum_ic->Initialize(lev, solid.momentum_mf, 0.0);
238 solid.energy_ic ->Initialize(lev, solid.energy_mf, 0.0);
239
240 ic_m0 ->Initialize(lev, m0_mf, 0.0);
241 ic_u0 ->Initialize(lev, u0_mf, 0.0);
242 ic_q ->Initialize(lev, q_mf, 0.0);
243
244 Source_mf[lev] ->setVal(0.0);
245
246 if (managed) { if (lev >= (int)mixed.size()) mixed.push_back(false);}
247 else Mix(lev);
248}
249
250void Hydro::Mix(int lev)
251{
252 if (managed && mixed[lev]) return;
253
254 for (amrex::MFIter mfi(*velocity_mf[lev], true); mfi.isValid(); ++mfi)
255 {
256 const amrex::Box& bx = mfi.growntilebox();
257
258 Set::Patch<const Set::Scalar> eta_patch = eta_old_mf->Patch(lev,mfi);
259
268 Set::Patch<const Set::Scalar> rho_solid = solid.density_mf.Patch(lev,mfi);
269 Set::Patch<const Set::Scalar> M_solid = solid.momentum_mf.Patch(lev,mfi);
270 Set::Patch<const Set::Scalar> E_solid = solid.energy_mf.Patch(lev,mfi);
274
275
276 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
277 {
278 Set::Scalar eta = invert ? 1.0-eta_patch(i,j,k)*eta_patch(i,j,k) : eta_patch(i,j,k);
279
280 // Initially compute primitives (T,P,u) from given initial conditions
281 // But from then on, compute them from mixed values to avoid zero T conditions
282 // Except velocity - keep velocity from fluid values only
283 gas.ComputeLocalFractions(rho, Y, X, i,j,k); // Get local mole/mass fractions from fluid densities
284 Set::Scalar density = gas.ComputeD(rho, i, j, k); // If a gas mixture, this will compute the mixture density
285 T(i,j,k) = gas.ComputeT(p(i,j,k), density, X, i, j, k);
286 Set::Scalar E_fluid = gas.ComputeE(density, density*v(i,j,k,0), density*v(i,j,k,1), T(i,j,k), X, i, j, k);
287
288 // Mix
289 M(i, j, k, 0) = (rho(i, j, k)*v(i, j, k, 0))*eta + M_solid(i, j, k, 0)*(1.0-eta);
290 M(i, j, k, 1) = (rho(i, j, k)*v(i, j, k, 1))*eta + M_solid(i, j, k, 1)*(1.0-eta);
291 M_old(i, j, k, 0) = M(i, j, k, 0);
292 M_old(i, j, k, 1) = M(i, j, k, 1);
293
294 rho(i, j, k) = eta * rho(i, j, k) + (1.0 - eta) * rho_solid(i, j, k);
295 rho_old(i, j, k) = rho(i, j, k);
296
297 E(i, j, k) = E_fluid*eta + E_solid(i,j,k)*(1.0-eta);
298 E_old(i, j, k) = E(i, j, k);
299 //Util::Message(INFO,"Energy: ", E(i,j,k), " Pressure: ", p(i,j,k), " Temp: ", T(i,j,k), " Density: ",density, " R: ", gas.R(X,i,j,k), " MW: ", gas.GetMW(X,i,j,k), " Rg: ", Set::Constant::Rg);
300
301 //gas.ComputeLocalFractions(rho, Y, X, i,j,k); // Get local mole/mass fractions from mixed densities
302 //density = gas.ComputeD(rho, i, j, k);
303 //T(i, j, k) = gas.ComputeT(density, M(i,j,k,0), M(i,j,k,1), E(i,j,k), T(i,j,k), X, i, j, k);
304 //p(i, j, k) = gas.ComputeP(density, T(i,j,k), X, i, j, k);
305 //v(i,j,k,0) = M(i,j,k,0)/density;
306 //v(i,j,k,1) = M(i,j,k,1)/density;
307 });
308 //Util::Abort(INFO);
309 }
310 c_max = 0.0;
311 vx_max = 0.0;
312 vy_max = 0.0;
313}
314
315void Hydro::UpdateEta(int lev, Set::Scalar time)
316{
317 Util::Assert(INFO,TEST(!managed),"Should override this if Hydro is managed!");
318 eta_ic->Initialize(lev, *eta_mf, time);
319}
320
321void Hydro::UpdateFluxes(int /*lev*/, Set::Scalar /*time*/, Set::Scalar /*dt*/)
322{
323 Util::Assert(INFO,TEST(!managed),"Should override this if Hydro is managed!");
324}
325
327{
328
329}
330
332{
333 if (dynamictimestep.on)
335 return;
336
337 const Set::Scalar* DX = geom[lev].CellSize();
338
339 amrex::ParallelDescriptor::ReduceRealMax(c_max);
340 amrex::ParallelDescriptor::ReduceRealMax(vx_max);
341 amrex::ParallelDescriptor::ReduceRealMax(vy_max);
342
343 Set::Scalar new_timestep = cfl / ((c_max + vx_max) / DX[0] + (c_max + vy_max) / DX[1]);
344
345 Util::Assert(INFO, TEST(AMREX_SPACEDIM == 2));
346
347 SetTimestep(new_timestep);
348}
349
351{
352
353 if (!managed) std::swap(*eta_old_mf, *eta_mf);
354 std::swap(density_old_mf[lev], density_mf[lev]);
355 std::swap(momentum_old_mf[lev], momentum_mf[lev]);
356 std::swap(energy_old_mf[lev], energy_mf[lev]);
357
358 //
359 // UPDATE ETA AND CALCULATE ETADOT
360 //
361
362 if (!managed) UpdateEta(lev, time);
363 if (managed)
364 {
365 UpdateFluxes(lev,time,dt);
366 Mix(lev);
367 }
368 for (amrex::MFIter mfi(*(velocity_mf)[lev], true); mfi.isValid(); ++mfi)
369 {
370 const amrex::Box& bx = mfi.growntilebox();
371 amrex::Array4<const Set::Scalar> const& eta_new = (*(*eta_mf)[lev]).array(mfi);
372 amrex::Array4<const Set::Scalar> const& eta = (*(*eta_old_mf)[lev]).array(mfi);
373 amrex::Array4<Set::Scalar> const& etadot = (*etadot_mf[lev]).array(mfi);
374 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
375 {
376
377 etadot(i, j, k) = (eta_new(i, j, k) - eta(i, j, k)) / dt;
378 if (invert) etadot(i,j,k) *= 1.0;
379
380 });
381 }
382
383
384 //
385 // DO TIME INTEGRATION (driving the RHS function)
386 //
387
388 // Organize references to the "new" solution
389 amrex::Vector<amrex::MultiFab> solution_new;
390 solution_new.emplace_back(*density_mf[lev].get(),amrex::MakeType::make_alias,0,1);
391 solution_new.emplace_back(*momentum_mf[lev].get(),amrex::MakeType::make_alias,0,2);
392 solution_new.emplace_back(*energy_mf[lev].get(),amrex::MakeType::make_alias,0,1);
393
394 // Organize references to the "old" solution
395 amrex::Vector<amrex::MultiFab> solution_old;
396 solution_old.emplace_back(*density_old_mf[lev].get(),amrex::MakeType::make_alias,0,1);
397 solution_old.emplace_back(*momentum_old_mf[lev].get(),amrex::MakeType::make_alias,0,2);
398 solution_old.emplace_back(*energy_old_mf[lev].get(),amrex::MakeType::make_alias,0,1);
399
400 // Create the time integrator
401 amrex::TimeIntegrator timeintegrator(solution_new, time);
402
403 // Set the time integrator RHS - in this case, just relay to our current RHS function
404 timeintegrator.set_rhs([&](amrex::Vector<amrex::MultiFab> & rhs_mf, amrex::Vector<amrex::MultiFab> & solution_mf, const Set::Scalar time)
405 {
406 RHS(lev, time,
407 rhs_mf[0], rhs_mf[1], rhs_mf[2],
408 solution_mf[0],solution_mf[1],solution_mf[2]);
409 });
410
411 // Take care of filling boundaries during stages
412 timeintegrator.set_post_stage_action([&](amrex::Vector<amrex::MultiFab> & stage_mf, Set::Scalar time)
413 {
414 density_bc->FillBoundary(stage_mf[0],0,1,time,0);
415 stage_mf[0].FillBoundary(true);
416 momentum_bc->FillBoundary(stage_mf[1],0,2,time,0);
417 stage_mf[1].FillBoundary(true);
418 energy_bc->FillBoundary(stage_mf[2],0,1,time,0);
419 stage_mf[2].FillBoundary(true);
420 });
421
422 // Do the update
423 timeintegrator.advance(solution_old, solution_new, time, dt);
424
425
426 //
427 // APPLY CUTOFFS AND DO DYNAMIC TIMESTEP CALCULATION
428 //
429
430 Set::Scalar dt_max = std::numeric_limits<Set::Scalar>::max();
431 for (amrex::MFIter mfi(*velocity_mf[lev], false); mfi.isValid(); ++mfi)
432 {
433 const amrex::Box& bx = mfi.validbox();
434 const Set::Scalar* DX = geom[lev].CellSize();
435
436 Set::Patch<const Set::Scalar> eta_patch = eta_mf->Patch(lev,mfi);
437 Set::Patch<const Set::Scalar> rho_solid = solid.density_mf.Patch(lev,mfi);
438 Set::Patch<const Set::Scalar> M_solid = solid.momentum_mf.Patch(lev,mfi);
439 Set::Patch<const Set::Scalar> E_solid = solid.energy_mf.Patch(lev,mfi);
440
441 Set::Patch<Set::Scalar> rho_new = density_mf.Patch(lev,mfi);
442 Set::Patch<Set::Scalar> E_new = energy_mf.Patch(lev,mfi);
444
446
448 Set::Patch<Set::Scalar> Source = Source_mf.Patch(lev,mfi);
449
450 Set::Scalar *dt_max_handle = &dt_max;
451
452 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
453 {
454 Set::Scalar eta = invert ? 1.0-eta_patch(i,j,k)*eta_patch(i,j,k) : eta_patch(i,j,k);
455
456 if (eta < cutoff)
457 {
458 rho_new(i,j,k,0) = rho_solid(i,j,k,0);
459 M_new(i,j,k,0) = M_solid(i,j,k,0);
460 M_new(i,j,k,1) = M_solid(i,j,k,1);
461 E_new(i,j,k,0) = E_solid(i,j,k,0);
462 }
463
464 Set::Matrix gradu = Numeric::Gradient(u, i, j, k, DX);
465 omega(i, j, k) = eta * (gradu(1,0) - gradu(0,1));
466
467 if (dynamictimestep.on)
468 {
469 *dt_max_handle = std::fabs(cfl * DX[0] / (u(i,j,k,0)*eta + small));
470 *dt_max_handle = std::min(*dt_max_handle, std::fabs(cfl * DX[1] / (u(i,j,k,1)*eta + small)));
471 *dt_max_handle = std::min(*dt_max_handle, std::fabs(cfl_v * DX[0]*DX[0] / (Source(i,j,k,1)+small)));
472 *dt_max_handle = std::min(*dt_max_handle, std::fabs(cfl_v * DX[1]*DX[1] / (Source(i,j,k,2)+small)));
473 }
474 });
475 }
476
477
478 if (dynamictimestep.on)
479 {
480 this->DynamicTimestep_SyncTimeStep(lev,dt_max);
481 }
482
483}//end Advance
484
485
486void Hydro::RHS(int lev, Set::Scalar /*time*/,
487 amrex::MultiFab &rho_rhs_mf,
488 amrex::MultiFab &M_rhs_mf,
489 amrex::MultiFab &E_rhs_mf,
490 const amrex::MultiFab &rho_mf,
491 const amrex::MultiFab &M_mf,
492 const amrex::MultiFab &E_mf)
493{
494
495 for (amrex::MFIter mfi(*(velocity_mf)[lev], true); mfi.isValid(); ++mfi)
496 {
497 const amrex::Box& bx = mfi.growntilebox();
498 amrex::Array4<const Set::Scalar> const& eta_patch = (*(*eta_old_mf)[lev]).array(mfi);
499
500 Set::Patch<const Set::Scalar> rho = rho_mf.array(mfi); // density
501 Set::Patch<const Set::Scalar> M = M_mf.array(mfi); // momentum
502 Set::Patch<const Set::Scalar> E = E_mf.array(mfi); // total energy (internal energy + kinetic energy) per unit volume (E/rho = e + 0.5*v^2)
503
504 Set::Patch<const Set::Scalar> rho_solid = solid.density_mf.Patch(lev,mfi);
505 Set::Patch<const Set::Scalar> M_solid = solid.momentum_mf.Patch(lev,mfi);
506
507 Set::Patch<Set::Scalar> scratch = scratch_mf.Patch(lev,mfi);
508
514
515 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
516 {
517 Set::Scalar eta = invert ? 1.0-eta_patch(i,j,k)*eta_patch(i,j,k) : eta_patch(i,j,k);
518
519 // Compute T and P primitives from mixed values
520 Set::Scalar density = gas.ComputeD(rho, i, j, k);
521 T(i,j,k) = gas.ComputeT(density, M(i,j,k,0), M(i,j,k,1), E(i,j,k), T(i,j,k), X, i, j, k);
522 p(i,j,k) = gas.ComputeP(density, T(i,j,k), X, i, j, k);
523
524 // Compute velocity from fluid values
525 scratch(i,j,k) = (rho(i,j,k) - rho_solid(i,j,k)*(1.0 - eta))/(eta + small);
526 gas.ComputeLocalFractions(scratch, Y, X, i, j, k);
527 Set::Scalar density_fluid = gas.ComputeD(scratch, i, j, k);
528 Set::Scalar Mx_fluid = (M(i,j,k,0) - M_solid(i,j,k,0)*(1.0 - eta))/(eta + small);
529 Set::Scalar My_fluid = (M(i,j,k,1) - M_solid(i,j,k,1)*(1.0 - eta))/(eta + small);
530 v(i,j,k,0) = Mx_fluid/density_fluid;
531 v(i,j,k,1) = My_fluid/density_fluid;
532
533 if (eta < small)
534 {
535 v(i,j,k,0) *= eta;
536 v(i,j,k,1) *= eta;
537
538 #if AMREX_SPACEDIM == 3
539 v(i,j,k,2) *= eta;
540 #endif
541 }
542 });
543 }
544
545 const Set::Scalar* DX = geom[lev].CellSize();
546 amrex::Box domain = geom[lev].Domain();
547
548 for (amrex::MFIter mfi(*(*eta_mf)[lev], false); mfi.isValid(); ++mfi)
549 {
550 const amrex::Box& bx = mfi.validbox();
551
552 // Inputs
553 Set::Patch<const Set::Scalar> rho = rho_mf.array(mfi);
554 Set::Patch<const Set::Scalar> E = E_mf.array(mfi);
555 Set::Patch<const Set::Scalar> M = M_mf.array(mfi);
556
557 // Outputs
558 Set::Patch<Set::Scalar> rho_rhs = rho_rhs_mf.array(mfi);
559 Set::Patch<Set::Scalar> M_rhs = M_rhs_mf.array(mfi);
560 Set::Patch<Set::Scalar> E_rhs = E_rhs_mf.array(mfi);
561
562
563 // Set::Patch<Set::Scalar> rho_new = density_mf.Patch(lev,mfi);
564 // Set::Patch<Set::Scalar> E_new = energy_mf.Patch(lev,mfi);
565 // Set::Patch<Set::Scalar> M_new = momentum_mf.Patch(lev,mfi);
566
567 Set::Patch<const Set::Scalar> rho_solid = solid.density_mf.Patch(lev,mfi);
568 Set::Patch<const Set::Scalar> M_solid = solid.momentum_mf.Patch(lev,mfi);
569 Set::Patch<const Set::Scalar> E_solid = solid.energy_mf.Patch(lev,mfi);
570
572
573 Set::Patch<const Set::Scalar> eta_patch = eta_old_mf->Patch(lev,mfi);
578
582
583 amrex::Array4<Set::Scalar> const& Source = (*Source_mf[lev]).array(mfi);
584
585 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
586 {
587 auto sten = Numeric::GetStencil(i, j, k, domain);
588
589 Set::Scalar eta = invert ? 1.0-eta_patch(i,j,k)*eta_patch(i,j,k) : eta_patch(i,j,k);
590
591 //Diffuse Sources
592 Set::Vector grad_eta = Numeric::Gradient(eta_patch, i, j, k, 0, DX);
593 Set::Scalar grad_eta_mag = grad_eta.lpNorm<2>();
594 Set::Matrix hess_eta = Numeric::Hessian(eta_patch, i, j, k, 0, DX);
595 if (invert) grad_eta *= -1.0;
596 if (invert) hess_eta *= -1.0;
597
598 #if AMREX_SPACEDIM == 2
599 Set::Vector u = Set::Vector(velocity(i, j, k, 0), velocity(i, j, k, 1)); // Velocity
600 Set::Vector u0 = Set::Vector(_u0(i, j, k, 0), _u0(i, j, k, 1)); // Velocity
601 Set::Vector q0 = Set::Vector(q(i,j,k,0), q(i,j,k,1));
602 #endif
603
604 #if AMREX_SPACEDIM == 3
605 Set::Vector u = Set::Vector(velocity(i, j, k, 0), velocity(i, j, k, 1), velocity(i, j, k, 2)); // Velocity
606 Set::Vector u0 = Set::Vector(_u0(i, j, k, 0), _u0(i, j, k, 1), _u0(i, j, k, 2)); // Velocity
607 Set::Vector q0 = Set::Vector(q(i,j,k,0), q(i,j,k,1), q(i,j,k,2));
608 #endif
609
610 Set::Matrix gradM = Numeric::Gradient(M, i, j, k, DX);
611 Set::Vector gradrho = Numeric::Gradient(rho,i,j,k,0,DX);
612 Set::Matrix hess_rho = Numeric::Hessian(rho,i,j,k,0,DX,sten);
613 Set::Matrix gradu = (gradM - u*gradrho.transpose()) / rho(i,j,k);
614
616 {
617 Set::Vector N = grad_eta / (grad_eta_mag + small);
618 // Set::Vector T(N(1), -N(0));
619 // u0 = N * u0(0) + T * u0(1);
620
621 #if AMREX_SPACEDIM == 2
622 Set::Vector T(N(1), -N(0));
623 u0 = N * u0(0) + T * u0(1);
624 #endif
625
626 #if AMREX_SPACEDIM == 3
627 Set::Vector T;
628 T(0) = N(1);
629 T(1) = -N(0);
630 T(2) = 0;
631 u0 = N*u0(0) + T * u0(1);
632 // Might not be physcially accurate, need to find how to extend to 3 dimensions
633 #endif
634 }
635
636
637 Set::Scalar mdot0 = m0(i,j,k)*grad_eta_mag;
638 Set::Vector Pdot0 = Set::Vector::Zero(); // Linear momentum source term
639 Set::Scalar qdot0 = q0.dot(grad_eta);
640
641 Set::Scalar mu = gas.dynamic_viscosity(T(i,j,k), molef, i, j, k);
642
643 // sten is necessary here because sometimes corner ghost
644 // cells don't get filled
645 Set::Matrix3 hess_M = Numeric::Hessian(M,i,j,k,DX);
647 for (int p = 0; p < 2; p++)
648 for (int q = 0; q < 2; q++)
649 for (int r = 0; r < 2; r++)
650 {
651 hess_u(r,p,q) =
652 (hess_M(r,p,q) - gradu(r,q)*gradrho(p) - gradu(r,p)*gradrho(q) - u(r)*hess_rho(p,q))
653 / rho(i,j,k);
654 }
655
656 Set::Vector Ldot0 = Set::Vector::Zero();
657 Set::Vector div_tau = Set::Vector::Zero();
658 Set::Scalar lambda = 0.0; //-2.0/3.0*mu_eff;
659 for (int p = 0; p<2; p++)
660 for (int q = 0; q<2; q++)
661 for (int r = 0; r<2; r++)
662 for (int s = 0; s<2; s++)
663 {
664 Ldot0(p) += 0.25 * (mu * ((p==r && q==s) + (p==s && q==r)) + lambda * (p==q && r==s)) * (u(r) - u0(r)) * hess_eta(q, s);
665 div_tau(p) += 0.5 * (mu * ((p==r && q==s) + (p==s && q==r)) + lambda * (p==q && r==s)) * (hess_u(r,q,s) + hess_u(s,q,r));
666
667 }
668
669 Source(i,j, k, 0) = mdot0;
670 Source(i,j, k, 1) = Pdot0(0) - Ldot0(0);
671 Source(i,j, k, 2) = Pdot0(1) - Ldot0(1);
672 Source(i,j, k, 3) = qdot0;// - Ldot0(0)*v(i,j,k,0) - Ldot0(1)*v(i,j,k,1);
673
674 // Lagrange terms to enforce no-penetration
675 Source(i,j,k,1) -= lagrange*(u-u0).dot(grad_eta)*grad_eta(0);
676 Source(i,j,k,2) -= lagrange*(u-u0).dot(grad_eta)*grad_eta(1);
677
678 //Godunov flux
679 //states of total fields
680 const int X = 0, Y = 1;
681 Solver::Local::Riemann::State state_xlo(rho, M, E, i-1, j, k, X);
682 Solver::Local::Riemann::State state_x (rho, M, E, i , j, k, X);
683 Solver::Local::Riemann::State state_xhi(rho, M, E, i+1, j, k, X);
684
685 Solver::Local::Riemann::State state_ylo(rho, M, E, i, j-1, k, Y);
686 Solver::Local::Riemann::State state_y (rho, M, E, i, j , k, Y);
687 Solver::Local::Riemann::State state_yhi(rho, M, E, i, j+1, k, Y);
688
689 //states of solid fields
690 Solver::Local::Riemann::State state_xlo_solid(rho_solid, M_solid, E_solid, i-1, j, k, X);
691 Solver::Local::Riemann::State state_x_solid (rho_solid, M_solid, E_solid, i , j, k, X);
692 Solver::Local::Riemann::State state_xhi_solid(rho_solid, M_solid, E_solid, i+1, j, k, X);
693
694 Solver::Local::Riemann::State state_ylo_solid(rho_solid, M_solid, E_solid, i, j-1, k, Y);
695 Solver::Local::Riemann::State state_y_solid (rho_solid, M_solid, E_solid, i, j , k, Y);
696 Solver::Local::Riemann::State state_yhi_solid(rho_solid, M_solid, E_solid, i, j+1, k, Y);
697
698 Solver::Local::Riemann::State state_xlo_fluid = invert ?
699 (state_xlo - (eta_patch(i-1,j,k))*state_xlo_solid) / (1.0 - eta_patch(i-1,j,k) + small) :
700 (state_xlo - (1.0 - eta_patch(i-1,j,k))*state_xlo_solid) / (eta_patch(i-1,j,k) + small);
701 Solver::Local::Riemann::State state_x_fluid = invert ?
702 (state_x - (eta_patch(i,j,k) )*state_x_solid ) / (1.0 - eta_patch(i,j,k) + small):
703 (state_x - (1.0 - eta_patch(i,j,k) )*state_x_solid ) / (eta_patch(i,j,k) + small);
704 Solver::Local::Riemann::State state_xhi_fluid = invert ?
705 (state_xhi - (eta_patch(i+1,j,k))*state_xhi_solid) / (1.0 - eta_patch(i+1,j,k) + small) :
706 (state_xhi - (1.0 - eta_patch(i+1,j,k))*state_xhi_solid) / (eta_patch(i+1,j,k) + small);
707 Solver::Local::Riemann::State state_ylo_fluid = invert ?
708 (state_ylo - (eta_patch(i,j-1,k))*state_ylo_solid) / (1.0 - eta_patch(i,j-1,k) + small):
709 (state_ylo - (1.0 - eta_patch(i,j-1,k))*state_ylo_solid) / (eta_patch(i,j-1,k) + small);
710 Solver::Local::Riemann::State state_y_fluid = invert ?
711 (state_y - (eta_patch(i,j,k) )*state_y_solid ) / (1.0 - eta_patch(i,j,k) + small):
712 (state_y - (1.0 - eta_patch(i,j,k) )*state_y_solid ) / (eta_patch(i,j,k) + small);
713 Solver::Local::Riemann::State state_yhi_fluid = invert ?
714 (state_yhi - (eta_patch(i,j+1,k))*state_yhi_solid) / (1.0 - eta_patch(i,j+1,k) + small):
715 (state_yhi - (1.0 - eta_patch(i,j+1,k))*state_yhi_solid) / (eta_patch(i,j+1,k) + small);
716
717 Solver::Local::Riemann::Flux flux_xlo, flux_ylo, flux_xhi, flux_yhi;
718
719 try
720 {
721 //lo interface fluxes
722 flux_xlo = riemannsolver->Solve(state_xlo_fluid, state_x_fluid, gas, molef, i, j, k, 0, small) * eta;
723 flux_ylo = riemannsolver->Solve(state_ylo_fluid, state_y_fluid, gas, molef, i, j, k, 2, small) * eta;
724
725 //hi interface fluxes
726 flux_xhi = riemannsolver->Solve(state_x_fluid, state_xhi_fluid, gas, molef, i, j, k, 1, small) * eta;
727 flux_yhi = riemannsolver->Solve(state_y_fluid, state_yhi_fluid, gas, molef, i, j, k, 3, small) * eta;
728 }
729 catch(...)
730 {
731 Util::ParallelMessage(INFO,"lev=",lev);
732 Util::ParallelMessage(INFO,"i=",i,"j=",j);
734 }
735
736
737 Set::Scalar drhof_dt =
738 (flux_xlo.mass - flux_xhi.mass) / DX[0] +
739 (flux_ylo.mass - flux_yhi.mass) / DX[1] +
740 Source(i, j, k, 0);
741
742 rho_rhs(i,j,k) =
743 // rho_new(i, j, k) = rho(i, j, k) +
744 //(
745 drhof_dt +
746 // todo add drhos_dt term if want time-evolving rhos
747 etadot(i,j,k) * (rho(i,j,k) - rho_solid(i,j,k)) / (eta + small)
748 // ) * dt;
749 ;
750
751
752
753 Set::Scalar dMxf_dt =
754 (flux_xlo.momentum_normal - flux_xhi.momentum_normal ) / DX[0] +
755 (flux_ylo.momentum_tangent - flux_yhi.momentum_tangent) / DX[1] +
756 div_tau(0) * eta +
757 g(0)*rho(i,j,k) +
758 Source(i, j, k, 1);
759
760 M_rhs(i,j,k,0) =
761 //M_new(i, j, k, 0) = M(i, j, k, 0) +
762 // (
763 dMxf_dt +
764 // todo add dMs_dt term if want time-evolving Ms
765 etadot(i,j,k)*(M(i,j,k,0) - M_solid(i,j,k,0)) / (eta + small)
766 // ) * dt;
767 ;
768
769 Set::Scalar dMyf_dt =
770 (flux_xlo.momentum_tangent - flux_xhi.momentum_tangent) / DX[0] +
771 (flux_ylo.momentum_normal - flux_yhi.momentum_normal ) / DX[1] +
772 div_tau(1) * eta +
773 g(1)*rho(i,j,k) +
774 Source(i, j, k, 2);
775
776 M_rhs(i,j,k,1) =
777 //M_new(i, j, k, 1) = M(i, j, k, 1) +
778 //(
779 dMyf_dt +
780 // todo add dMs_dt term if want time-evolving Ms
781 etadot(i,j,k)*(M(i,j,k,1) - M_solid(i,j,k,1)) / (eta+small)
782 // )*dt;
783 ;
784
785 Set::Scalar dEf_dt =
786 (flux_xlo.energy - flux_xhi.energy) / DX[0] +
787 (flux_ylo.energy - flux_yhi.energy) / DX[1] +
788 Source(i, j, k, 3);
789
790 E_rhs(i,j,k) =
791 // E_new(i, j, k) = E(i, j, k) +
792 // (
793 dEf_dt +
794 // todo add dEs_dt term if want time-evolving Es
795 etadot(i,j,k)*(E(i,j,k) - E_solid(i,j,k)) / (eta+small)
796 // ) * dt;
797 ;
798
799#ifdef AMREX_DEBUG
800 if ((rho_rhs(i,j,k) != rho_rhs(i,j,k)) ||
801 (M_rhs(i,j,k,0) != M_rhs(i,j,k,0)) ||
802 (M_rhs(i,j,k,1) != M_rhs(i,j,k,1)) ||
803 (E_rhs(i,j,k) != E_rhs(i,j,k)))
804 {
805 Util::ParallelMessage(INFO,"rho_rhs=",rho_rhs(i,j,k));
806 Util::ParallelMessage(INFO,"Mx_rhs=",M_rhs(i,j,k,0));
807 Util::ParallelMessage(INFO,"Mx_rhs=",M_rhs(i,j,k,1));
808 Util::ParallelMessage(INFO,"E_rhs=",E_rhs(i,j,k));
809
810 Util::ParallelMessage(INFO,"lev=",lev);
811 Util::ParallelMessage(INFO,"i=",i," j=",j);
812 Util::ParallelMessage(INFO,"drhof_dt ",drhof_dt); // dies
813 Util::ParallelMessage(INFO,"flux_xlo.mass ",flux_xlo.mass);
814 Util::ParallelMessage(INFO,"flux_xhi.mass ",flux_xhi.mass); // dies, depends on state_xx, state_xhi, state_x_solid, state_xhi_solid, eta, small
815 Util::ParallelMessage(INFO,"flux_ylo.mass ",flux_ylo.mass);
816 Util::ParallelMessage(INFO,"flux_xhi.mass ",flux_yhi.mass);
817 Util::ParallelMessage(INFO,"eta ",eta);
818 Util::ParallelMessage(INFO,"etadot ",etadot(i,j,k));
819 Util::ParallelMessage(INFO,"Source ",Source(i,j,k,0));
820 Util::ParallelMessage(INFO,"state_x ",state_x); // <<<<
821 Util::ParallelMessage(INFO,"state_y ",state_y);
822 Util::ParallelMessage(INFO,"state_x_solid ",state_x_solid); // <<<<
823 Util::ParallelMessage(INFO,"state_y_solid ",state_y_solid);
824 Util::ParallelMessage(INFO,"state_xhi ",state_xhi); // <<<<
825 Util::ParallelMessage(INFO,"state_yhi ",state_yhi);
826 Util::ParallelMessage(INFO,"state_xhi_solid ",state_xhi_solid);
827 Util::ParallelMessage(INFO,"state_yhi_solids ",state_yhi_solid);
828 Util::ParallelMessage(INFO,"state_xlo ",state_xlo);
829 Util::ParallelMessage(INFO,"state_ylo ",state_ylo);
830 Util::ParallelMessage(INFO,"state_xlo_solid ",state_xlo_solid);
831 Util::ParallelMessage(INFO,"state_ylo_solid ",state_ylo_solid);
832
833 Util::ParallelMessage(INFO,"Mx_solid ",M_solid(i,j,k,0));
834 Util::ParallelMessage(INFO,"My_solid ",M_solid(i,j,k,1));
836 Util::ParallelMessage(INFO,"Mx ",M(i,j,k,0));
837 Util::ParallelMessage(INFO,"My ",M(i,j,k,1));
838 Util::ParallelMessage(INFO,"dMx/dt ",dMxf_dt);
839 Util::ParallelMessage(INFO,"dMy/dt ",dMyf_dt);
840
841
844 Util::Message(INFO,DX[0]);
847 Util::Message(INFO,DX[1]);
848 Util::Message(INFO,div_tau);
849 Util::Message(INFO,Source(i, j, k, 2));
850
851 Util::Message(INFO,hess_eta);
852 Util::Message(INFO,velocity(i,j,k,0));
853 Util::Message(INFO,velocity(i,j,k,1));
854
856 }
857#endif
858
859
860
861 // todo - may need to move this for higher order schemes...
862 omega(i, j, k) = eta * (gradu(1,0) - gradu(0,1));
863 });
864 }
865}
866
867void Hydro::Regrid(int lev, Set::Scalar /* time */)
868{
869 BL_PROFILE("Integrator::Hydro::Regrid");
870 Source_mf[lev]->setVal(0.0);
871 if (lev < finest_level) return;
872
873 Util::Message(INFO, "Regridding on level", lev);
874}//end regrid
875
876//void Hydro::TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar time, int ngrow)
877void Hydro::TagCellsForRefinement(int lev, amrex::TagBoxArray& a_tags, Set::Scalar, int)
878{
879 BL_PROFILE("Integrator::Flame::TagCellsForRefinement");
880
881 const Set::Scalar* DX = geom[lev].CellSize();
882 Set::Scalar dr = sqrt(AMREX_D_TERM(DX[0] * DX[0], +DX[1] * DX[1], +DX[2] * DX[2]));
883
884 // Eta criterion for refinement
885 for (amrex::MFIter mfi(*(*eta_mf)[lev], true); mfi.isValid(); ++mfi) {
886 const amrex::Box& bx = mfi.tilebox();
887 amrex::Array4<char> const& tags = a_tags.array(mfi);
888 amrex::Array4<const Set::Scalar> const& eta = (*(*eta_mf)[lev]).array(mfi);
889
890 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
891 Set::Vector grad_eta = Numeric::Gradient(eta, i, j, k, 0, DX);
892 if (grad_eta.lpNorm<2>() * dr * 2 > eta_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
893 });
894 }
895
896 // Vorticity criterion for refinement
897 for (amrex::MFIter mfi(*vorticity_mf[lev], true); mfi.isValid(); ++mfi) {
898 const amrex::Box& bx = mfi.tilebox();
899 amrex::Array4<char> const& tags = a_tags.array(mfi);
900 amrex::Array4<const Set::Scalar> const& omega = (*vorticity_mf[lev]).array(mfi);
901
902 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
903 auto sten = Numeric::GetStencil(i, j, k, bx);
904 Set::Vector grad_omega = Numeric::Gradient(omega, i, j, k, 0, DX, sten);
905 if (grad_omega.lpNorm<2>() * dr * 2 > omega_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
906 });
907 }
908
909 // Gradu criterion for refinement
910 for (amrex::MFIter mfi(*velocity_mf[lev], true); mfi.isValid(); ++mfi) {
911 const amrex::Box& bx = mfi.tilebox();
912 amrex::Array4<char> const& tags = a_tags.array(mfi);
913 amrex::Array4<const Set::Scalar> const& v = (*velocity_mf[lev]).array(mfi);
914
915 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
916 auto sten = Numeric::GetStencil(i, j, k, bx);
917 Set::Matrix grad_u = Numeric::Gradient(v, i, j, k, DX, sten);
918 if (grad_u.lpNorm<2>() * dr * 2 > gradu_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
919 });
920 }
921
922 // Pressure criterion for refinement
923 for (amrex::MFIter mfi(*pressure_mf[lev], true); mfi.isValid(); ++mfi) {
924 const amrex::Box& bx = mfi.tilebox();
925 amrex::Array4<char> const& tags = a_tags.array(mfi);
926 amrex::Array4<const Set::Scalar> const& p = (*pressure_mf[lev]).array(mfi);
927
928 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
929 auto sten = Numeric::GetStencil(i, j, k, bx);
930 Set::Vector grad_p = Numeric::Gradient(p, i, j, k, 0, DX, sten);
931 if (grad_p.lpNorm<2>() * dr * 2 > p_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
932 });
933 }
934
935 // Density criterion for refinement
936 for (amrex::MFIter mfi(*density_mf[lev], true); mfi.isValid(); ++mfi) {
937 const amrex::Box& bx = mfi.tilebox();
938 amrex::Array4<char> const& tags = a_tags.array(mfi);
939 amrex::Array4<const Set::Scalar> const& rho = (*density_mf[lev]).array(mfi);
940
941 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
942 auto sten = Numeric::GetStencil(i, j, k, bx);
943 Set::Vector grad_rho = Numeric::Gradient(rho, i, j, k, 0, DX, sten);
944 if (grad_rho.lpNorm<2>() * dr * 2 > rho_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
945 });
946 }
947
948}
949
950}
#define X(name)
#define pp_query_required(...)
Definition ParmParse.H:101
#define pp_query_default(...)
Definition ParmParse.H:102
#define pp_forbid(...)
Definition ParmParse.H:110
#define pp_queryclass(...)
Definition ParmParse.H:109
#define TEST(x)
Definition Util.H:25
#define INFO
Definition Util.H:24
virtual void FillBoundary(amrex::BaseFab< T > &in, const amrex::Box &box, int ngrow, int dcomp, int ncomp, amrex::Real time, Orientation face=Orientation::All, const amrex::Mask *mask=nullptr)=0
Definition BMP.H:22
void Initialize(const int &a_lev, Set::Field< T > &a_field, Set::Scalar a_time=0.0)
Definition IC.H:39
Initialize Laminates in a matrix.
Definition Laminate.H:16
Definition PNG.H:26
int queryarr_default(std::string name, std::vector< std::string > &value, std::vector< std::string > defaultvalue)
Definition ParmParse.H:660
void forbid(std::string name, std::string explanation)
Definition ParmParse.H:160
int AnyUnusedInputs(bool inscopeonly=true, bool verbose=false)
Definition ParmParse.H:895
void queryclass(std::string name, T *value)
Definition ParmParse.H:960
void select_default(std::string name, PTRTYPE *&ic_eta, Args &&... args)
Definition ParmParse.H:1095
static int AllUnusedInputs()
Definition ParmParse.H:933
int query_default(std::string name, T &value, T defaultvalue)
Definition ParmParse.H:293
int query_validate(std::string name, int &value, std::vector< int > possibleintvals)
Definition ParmParse.H:336
void Regrid(int lev, Set::Scalar time) override
Definition Hydro.cpp:867
Set::Scalar vy_max
Definition Hydro.H:152
void Mix(int lev)
Definition Hydro.cpp:250
IC::IC< Set::Scalar > * momentum_ic
Definition Hydro.H:111
Set::Scalar eta_refinement_criterion
Definition Hydro.H:154
std::vector< bool > mixed
Definition Hydro.H:186
Set::Field< Set::Scalar > etadot_mf
Definition Hydro.H:117
IC::IC< Set::Scalar > * ic_q
Definition Hydro.H:140
Model::Gas::Gas gas
Definition Hydro.H:167
Set::Field< Set::Scalar > scratch_mf
Definition Hydro.H:104
IC::IC< Set::Scalar > * ic_u0
Definition Hydro.H:139
PrescribedFlowMode prescribedflowmode
Definition Hydro.H:178
Set::Scalar gradu_refinement_criterion
Definition Hydro.H:154
Set::Scalar cfl_v
Definition Hydro.H:155
struct Integrator::Hydro::@7 solid
Set::Field< Set::Scalar > vorticity_mf
Definition Hydro.H:119
Set::Field< Set::Scalar > Source_mf
Definition Hydro.H:124
Set::Scalar cutoff
Definition Hydro.H:158
void TimeStepBegin(Set::Scalar a_time, int a_iter) override
Definition Hydro.cpp:326
Set::Scalar rho_refinement_criterion
Definition Hydro.H:154
IC::IC< Set::Scalar > * ic_m0
Definition Hydro.H:138
Set::Field< Set::Scalar > * eta_old_mf
Definition Hydro.H:116
void Initialize(int lev) override
Definition Hydro.cpp:217
Set::Field< Set::Scalar > momentum_mf
Definition Hydro.H:95
Set::Field< Set::Scalar > mole_fraction_mf
Definition Hydro.H:103
Set::Field< Set::Scalar > energy_old_mf
Definition Hydro.H:93
void TimeStepComplete(Set::Scalar time, int lev) override
Definition Hydro.cpp:331
void TagCellsForRefinement(int lev, amrex::TagBoxArray &tags, amrex::Real, int) override
Definition Hydro.cpp:877
Set::Field< Set::Scalar > velocity_mf
Definition Hydro.H:98
Set::Scalar c_max
Definition Hydro.H:150
virtual void UpdateFluxes(int lev, Set::Scalar time, Set::Scalar dt)
Definition Hydro.cpp:321
Set::Scalar cfl
Definition Hydro.H:155
BC::Constant neumann_bc_D
Definition Hydro.H:148
Set::Field< Set::Scalar > density_old_mf
Definition Hydro.H:90
BC::Constant neumann_bc_1
Definition Hydro.H:147
static void Parse(Hydro &value, IO::ParmParse &pp)
Definition Hydro.cpp:35
BC::Nothing bc_nothing
Definition Hydro.H:145
void RHS(int lev, Set::Scalar time, amrex::MultiFab &rho_rhs_mf, amrex::MultiFab &M_rhs_mf, amrex::MultiFab &E_rhs_mf, const amrex::MultiFab &rho_mf, const amrex::MultiFab &M_mf, const amrex::MultiFab &E_mf)
Definition Hydro.cpp:486
IC::IC< Set::Scalar > * pressure_ic
Definition Hydro.H:134
BC::BC< Set::Scalar > * density_bc
Definition Hydro.H:127
Set::Field< Set::Scalar > temperature_mf
Definition Hydro.H:100
Set::Scalar vx_max
Definition Hydro.H:151
IC::IC< Set::Scalar > * velocity_ic
Definition Hydro.H:133
Set::Scalar p_refinement_criterion
Definition Hydro.H:154
BC::BC< Set::Scalar > * eta_bc
Definition Hydro.H:130
Set::Field< Set::Scalar > q_mf
Definition Hydro.H:123
IC::IC< Set::Scalar > * density_ic
Definition Hydro.H:110
BC::BC< Set::Scalar > * energy_bc
Definition Hydro.H:129
Set::Scalar small
Definition Hydro.H:157
Set::Field< Set::Scalar > u0_mf
Definition Hydro.H:122
Set::Vector g
Definition Hydro.H:181
Set::Field< Set::Scalar > momentum_old_mf
Definition Hydro.H:96
Set::Scalar omega_refinement_criterion
Definition Hydro.H:154
BC::BC< Set::Scalar > * momentum_bc
Definition Hydro.H:128
Set::Field< Set::Scalar > m0_mf
Definition Hydro.H:121
Set::Scalar lagrange
Definition Hydro.H:159
Set::Field< Set::Scalar > * eta_mf
Definition Hydro.H:115
Set::Field< Set::Scalar > mass_fraction_mf
Definition Hydro.H:102
virtual void UpdateEta(int lev, Set::Scalar time)
Definition Hydro.cpp:315
Set::Field< Set::Scalar > energy_mf
Definition Hydro.H:92
Set::Field< Set::Scalar > density_mf
Definition Hydro.H:89
void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
Definition Hydro.cpp:350
IC::IC< Set::Scalar > * energy_ic
Definition Hydro.H:112
Set::Field< Set::Scalar > pressure_mf
Definition Hydro.H:99
Solver::Local::Riemann::Riemann * riemannsolver
Definition Hydro.H:162
IC::IC< Set::Scalar > * eta_ic
Definition Hydro.H:142
void DynamicTimestep_SyncTimeStep(int lev, Set::Scalar dt_min)
Params for the dynamic timestp.
Definition Integrator.H:271
void RegisterNewFab(Set::Field< Set::Scalar > &new_fab, BC::BC< Set::Scalar > *new_bc, int ncomp, int nghost, std::string name, bool writeout, bool evolving=true, std::vector< std::string > suffix={})
Add a new cell-based scalar field.
struct Integrator::Integrator::@8 dynamictimestep
amrex::Vector< amrex::Real > dt
Timesteps for each level of refinement.
Definition Integrator.H:394
void SetTimestep(Set::Scalar _timestep)
Utility to set the coarse-grid timestep.
virtual const char * model_name() const =0
double ComputeP(double density, double T, Set::Patch< const Set::Scalar > &X, int i, int j, int k) const
Definition Gas.cpp:66
double ComputeT(double density, double momentumx, double momentumy, double E, double Tguess, Set::Patch< const Set::Scalar > &X, int i, int j, int k, double rtol=1e-12) const
Definition Gas.cpp:52
Thermo::Thermo * thermo
Definition Gas.H:23
void ComputeLocalFractions(Set::Patch< const Set::Scalar > &density_mf, Set::Patch< Set::Scalar > &mass_fraction_mf, Set::Patch< Set::Scalar > &mole_fraction_mf, const int i, const int j, const int k)
Definition Gas.H:81
double dynamic_viscosity(double T, Set::Patch< const Set::Scalar > &X, int i, int j, int k) const
Definition Gas.cpp:35
double ComputeD(Set::Patch< const Set::Scalar > &rhoY, int i, int j, int k) const
Definition Gas.H:107
EOS::EOS * eos
Definition Gas.H:25
double ComputeE(double density, double momentumx, double momentumy, double T, Set::Patch< const Set::Scalar > &X, int i, int j, int k) const
Definition Gas.cpp:71
Transport::Transport * transport
Definition Gas.H:24
int nspecies
Definition Gas.H:21
virtual const char * model_name() const =0
virtual const char * model_name() const =0
amrex::Array4< Set::Scalar > Patch(const int lev, const amrex::MFIter &mfi) const &
Definition Set.H:263
static Matrix3 Zero()
Definition Matrix3.H:42
virtual Flux Solve(State lo, State hi, Model::Gas::Gas &gas, Set::Patch< const Set::Scalar > &X, int i, int j, int k, int side, Set::Scalar small)=0
Roe Riemann Solver based on Gas Dynamics - Culbert B. Laney.
Definition Roe.H:30
Collection of numerical integrator objects.
Definition AllenCahn.H:42
static AMREX_FORCE_INLINE std::array< StencilType, AMREX_SPACEDIM > GetStencil(const int i, const int j, const int k, const amrex::Box domain)
Definition Stencil.H:45
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)
Definition Stencil.H:1043
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)
Definition Stencil.H:681
amrex::Real Scalar
Definition Base.H:18
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition Base.H:19
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
Definition Base.H:22
void ParallelMessage(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:185
void Abort(const char *msg)
Definition Util.cpp:268
AMREX_FORCE_INLINE void Assert(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)
Definition Util.H:58
void Warning(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:202
void Message(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:129
void Exception(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:226
Set::Scalar momentum_normal
Definition Riemann.H:94
Set::Scalar momentum_tangent
Definition Riemann.H:95