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