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
17#include "AMReX_TimeIntegrator.H"
18
19#if AMREX_SPACEDIM == 2
20
21namespace Integrator
22{
23
24Hydro::Hydro(IO::ParmParse& pp) : Hydro()
25{
26 pp_queryclass(*this);
27}
28
29void
30Hydro::Parse(Hydro& value, IO::ParmParse& pp)
31{
32 BL_PROFILE("Integrator::Hydro::Hydro()");
33 {
34 // pp.query_default("r_refinement_criterion", value.r_refinement_criterion , 0.01);
35 // energy-based refinement
36 // pp.query_default("e_refinement_criterion", value.e_refinement_criterion , 0.01);
37 // momentum-based refinement
38 // pp.query_default("m_refinement_criterion", value.m_refinement_criterion , 0.01);
39
40 pp.forbid("scheme","use integration.type instead");
41
42 // eta-based refinement
43 pp.query_default("eta_refinement_criterion", value.eta_refinement_criterion , 0.01);
44 // vorticity-based refinement
45 pp.query_default("omega_refinement_criterion", value.omega_refinement_criterion, 0.01);
46 // velocity gradient-based refinement
47 pp.query_default("gradu_refinement_criterion", value.gradu_refinement_criterion, 0.01);
48 // pressure-based refinement
49 pp.query_default("p_refinement_criterion", value.p_refinement_criterion, 1e100);
50 // density-based refinement
51 pp.query_default("rho_refinement_criterion", value.rho_refinement_criterion, 1e100);
52
53 pp_query_required("gamma", value.gamma); // gamma for gamma law
54 pp_query_required("cfl", value.cfl); // cfl condition
55 pp_query_default("cfl_v", value.cfl_v,1E100); // cfl condition
56 pp_query_required("mu", value.mu); // linear viscosity coefficient
57 pp_forbid("Lfactor","replaced with mu");
58 //pp_query_default("Lfactor", value.Lfactor,1.0); // (to be removed) test factor for viscous source
59 pp_forbid("Pfactor","replaced with mu");
60 //pp_query_default("Pfactor", value.Pfactor,1.0); // (to be removed) test factor for viscous source
61 pp_query_default("pref", value.pref,1.0); // reference pressure for Roe solver
62
63 pp_forbid("rho.bc","--> density.bc");
64 pp_forbid("p.bc","--> pressure.bc");
65 pp_forbid("v.bc", "--> velocity.bc");
66 pp_forbid("pressure.bc","--> energy.bc");
67 pp_forbid("velocity.bc","--> momentum.bc");
68
69 // Boundary condition for density
70 pp.select_default<BC::Constant,BC::Expression>("density.bc",value.density_bc,1);
71 // Boundary condition for energy
72 pp.select_default<BC::Constant,BC::Expression>("energy.bc",value.energy_bc,1);
73 // Boundary condition for momentum
74 pp.select_default<BC::Constant,BC::Expression>("momentum.bc",value.momentum_bc,2);
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 pp_query_default("small",value.small,1E-8); // small regularization value
79 pp_query_default("cutoff",value.cutoff,-1E100); // cutoff value
80 pp_query_default("lagrange",value.lagrange,0.0); // lagrange no-penetration factor
81
82 pp_forbid("roefix","--> solver.roe.entropy_fix"); // Roe solver entropy fix
83
84 }
85 // Register FabFields:
86 {
87 int nghost = 1;
88
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);
92
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);
95
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);
98
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);
101
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);
105
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"});
109
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);
113
114 value.RegisterNewFab(value.Source_mf, &value.bc_nothing, 4, 0, "Source", true, false);
115 }
116
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");
125
126
127 // ORDER PARAMETER
128
129 // eta initial condition
130 pp.select_default<IC::Constant,IC::Laminate,IC::Expression,IC::BMP,IC::PNG>("eta.ic",value.eta_ic,value.geom);
131
132 // PRIMITIVE FIELD INITIAL CONDITIONS
133
134 // velocity initial condition
135 pp.select_default<IC::Constant,IC::Expression>("velocity.ic",value.velocity_ic,value.geom);
136 // solid pressure initial condition
137 pp.select_default<IC::Constant,IC::Expression>("pressure.ic",value.pressure_ic,value.geom);
138 // density initial condition type
139 pp.select_default<IC::Constant,IC::Expression>("density.ic",value.density_ic,value.geom);
140
141
142 // SOLID FIELDS
143
144 // solid momentum initial condition
145 pp.select_default<IC::Constant,IC::Expression>("solid.momentum.ic",value.solid.momentum_ic,value.geom);
146 // solid density initial condition
147 pp.select_default<IC::Constant,IC::Expression>("solid.density.ic",value.solid.density_ic,value.geom);
148 // solid energy initial condition
149 pp.select_default<IC::Constant,IC::Expression>("solid.energy.ic",value.solid.energy_ic,value.geom);
150
151
152 // DIFFUSE BOUNDARY SOURCES
153
154 // diffuse boundary prescribed mass flux
155 pp.select_default<IC::Constant,IC::Expression>("m0.ic",value.ic_m0,value.geom);
156 // diffuse boundary prescribed velocity
157 pp.select_default<IC::Constant,IC::Expression>("u0.ic",value.ic_u0,value.geom);
158 // diffuse boundary prescribed heat flux
159 pp.select_default<IC::Constant,IC::Expression>("q.ic",value.ic_q,value.geom);
160
161 // Riemann solver
164 Solver::Local::Riemann::HLLC>("solver",value.riemannsolver);
165
166
167
168 bool allow_unused;
169 // Set this to true to allow unused inputs without error.
170 // (Not recommended.)
171 pp.query_default("allow_unused",allow_unused,false);
172 if (!allow_unused && pp.AnyUnusedInputs(true, false))
173 {
174 Util::Warning(INFO,"The following inputs were specified but not used:");
175 pp.AllUnusedInputs();
176 Util::Exception(INFO,"Aborting. Specify 'allow_unused=True` to ignore this error.");
177 }
178}
179
180
181void Hydro::Initialize(int lev)
182{
183 BL_PROFILE("Integrator::Hydro::Initialize");
184
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);
188
189 //flux_mf[lev] ->setVal(0.0);
190
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);
194
195 density_ic ->Initialize(lev, density_old_mf, 0.0);
196
197
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);
201
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);
205
206 Source_mf[lev] ->setVal(0.0);
207
208 Mix(lev);
209}
210
211void Hydro::Mix(int lev)
212{
213 for (amrex::MFIter mfi(*eta_mf[lev], true); mfi.isValid(); ++mfi)
214 {
215 const amrex::Box& bx = mfi.growntilebox();
216
217 Set::Patch<const Set::Scalar> eta = eta_mf.Patch(lev,mfi);
218
219 Set::Patch<const Set::Scalar> v = velocity_mf.Patch(lev,mfi);
220 Set::Patch<const Set::Scalar> p = pressure_mf.Patch(lev,mfi);
221 Set::Patch<Set::Scalar> rho = density_mf.Patch(lev,mfi);
222 Set::Patch<Set::Scalar> rho_old = density_old_mf.Patch(lev,mfi);
223 Set::Patch<Set::Scalar> M = momentum_mf.Patch(lev,mfi);
224 Set::Patch<Set::Scalar> M_old = momentum_old_mf.Patch(lev,mfi);
225 Set::Patch<Set::Scalar> E = energy_mf.Patch(lev,mfi);
226 Set::Patch<Set::Scalar> E_old = energy_old_mf.Patch(lev,mfi);
227 Set::Patch<const Set::Scalar> rho_solid = solid.density_mf.Patch(lev,mfi);
228 Set::Patch<const Set::Scalar> M_solid = solid.momentum_mf.Patch(lev,mfi);
229 Set::Patch<const Set::Scalar> E_solid = solid.energy_mf.Patch(lev,mfi);
230
231
232 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
233 {
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);
236
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);
241
242 E(i, j, k) =
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)
244 +
245 E_solid(i, j, k) * (1.0 - eta(i, j, k));
246 E_old(i, j, k) = E(i, j, k);
247 });
248 }
249 c_max = 0.0;
250 vx_max = 0.0;
251 vy_max = 0.0;
252}
253
254void Hydro::UpdateEta(int lev, Set::Scalar time)
255{
256 eta_ic->Initialize(lev, eta_mf, time);
257}
258
259void Hydro::TimeStepBegin(Set::Scalar, int /*iter*/)
260{
261
262}
263
264void Hydro::TimeStepComplete(Set::Scalar, int lev)
265{
266 if (dynamictimestep.on)
267 Integrator::DynamicTimestep_Update();
268 return;
269
270 const Set::Scalar* DX = geom[lev].CellSize();
271
272 amrex::ParallelDescriptor::ReduceRealMax(c_max);
273 amrex::ParallelDescriptor::ReduceRealMax(vx_max);
274 amrex::ParallelDescriptor::ReduceRealMax(vy_max);
275
276 Set::Scalar new_timestep = cfl / ((c_max + vx_max) / DX[0] + (c_max + vy_max) / DX[1]);
277
278 Util::Assert(INFO, TEST(AMREX_SPACEDIM == 2));
279
280 SetTimestep(new_timestep);
281}
282
283void Hydro::Advance(int lev, Set::Scalar time, Set::Scalar dt)
284{
285
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]);
290
291 //
292 // UPDATE ETA AND CALCULATE ETADOT
293 //
294
295 UpdateEta(lev, time);
296 for (amrex::MFIter mfi(*eta_mf[lev], true); mfi.isValid(); ++mfi)
297 {
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)
303 {
304 etadot(i, j, k) = (eta_new(i, j, k) - eta(i, j, k)) / dt;
305 });
306 }
307
308
309 //
310 // DO TIME INTEGRATION (driving the RHS function)
311 //
312
313 // Organize references to the "new" solution
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);
318
319 // Organize references to the "old" solution
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);
324
325 // Create the time integrator
326 amrex::TimeIntegrator timeintegrator(solution_new, time);
327
328 // Set the time integrator RHS - in this case, just relay to our current RHS function
329 timeintegrator.set_rhs([&](amrex::Vector<amrex::MultiFab> & rhs_mf, amrex::Vector<amrex::MultiFab> & solution_mf, const Set::Scalar time)
330 {
331 RHS(lev, time,
332 rhs_mf[0], rhs_mf[1], rhs_mf[2],
333 solution_mf[0],solution_mf[1],solution_mf[2]);
334 });
335
336 // Take care of filling boundaries during stages
337 timeintegrator.set_post_stage_action([&](amrex::Vector<amrex::MultiFab> & stage_mf, Set::Scalar time)
338 {
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);
345 });
346
347 // Do the update
348 timeintegrator.advance(solution_old, solution_new, time, dt);
349
350
351 //
352 // APPLY CUTOFFS AND DO DYNAMIC TIMESTEP CALCULATION
353 //
354
355 Set::Scalar dt_max = std::numeric_limits<Set::Scalar>::max();
356 for (amrex::MFIter mfi(*eta_mf[lev], false); mfi.isValid(); ++mfi)
357 {
358 const amrex::Box& bx = mfi.validbox();
359 const Set::Scalar* DX = geom[lev].CellSize();
360
361 Set::Patch<const Set::Scalar> eta = eta_mf.Patch(lev,mfi);
362 Set::Patch<const Set::Scalar> rho_solid = solid.density_mf.Patch(lev,mfi);
363 Set::Patch<const Set::Scalar> M_solid = solid.momentum_mf.Patch(lev,mfi);
364 Set::Patch<const Set::Scalar> E_solid = solid.energy_mf.Patch(lev,mfi);
365
366 Set::Patch<Set::Scalar> rho_new = density_mf.Patch(lev,mfi);
367 Set::Patch<Set::Scalar> E_new = energy_mf.Patch(lev,mfi);
368 Set::Patch<Set::Scalar> M_new = momentum_mf.Patch(lev,mfi);
369
370 Set::Patch<Set::Scalar> omega = vorticity_mf.Patch(lev,mfi);
371
372 Set::Patch<Set::Scalar> u = velocity_mf.Patch(lev,mfi);
373 Set::Patch<Set::Scalar> Source = Source_mf.Patch(lev,mfi);
374
375 Set::Scalar *dt_max_handle = &dt_max;
376
377 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
378 {
379 if (eta(i,j,k) < cutoff)
380 {
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);
385 }
386
387 Set::Matrix gradu = Numeric::Gradient(u, i, j, k, DX);
388 omega(i, j, k) = eta(i, j, k) * (gradu(1,0) - gradu(0,1));
389
390 if (dynamictimestep.on)
391 {
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)));
396 }
397 });
398 }
399
400
401 if (dynamictimestep.on)
402 {
403 this->DynamicTimestep_SyncTimeStep(lev,dt_max);
404 }
405
406}//end Advance
407
408
409void Hydro::RHS(int lev, Set::Scalar /*time*/,
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)
416{
417
418 for (amrex::MFIter mfi(*eta_mf[lev], true); mfi.isValid(); ++mfi)
419 {
420 const amrex::Box& bx = mfi.growntilebox();
421 amrex::Array4<const Set::Scalar> const& eta = (*eta_old_mf[lev]).array(mfi);
422
423 Set::Patch<const Set::Scalar> rho = rho_mf.array(mfi); // density
424 Set::Patch<const Set::Scalar> M = M_mf.array(mfi); // momentum
425 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)
426
427 Set::Patch<const Set::Scalar> rho_solid = solid.density_mf.Patch(lev,mfi);
428 Set::Patch<const Set::Scalar> M_solid = solid.momentum_mf.Patch(lev,mfi);
429 Set::Patch<const Set::Scalar> E_solid = solid.energy_mf.Patch(lev,mfi);
430
431 Set::Patch<Set::Scalar> v = velocity_mf.Patch(lev,mfi);
432 Set::Patch<Set::Scalar> p = pressure_mf.Patch(lev,mfi);
433
434 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
435 {
436
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);
439
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) );
442
443 //THESE ARE FLUID VELOCITY AND PRESSURE
444
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;
448 });
449 }
450
451 const Set::Scalar* DX = geom[lev].CellSize();
452 amrex::Box domain = geom[lev].Domain();
453
454 for (amrex::MFIter mfi(*eta_mf[lev], false); mfi.isValid(); ++mfi)
455 {
456 const amrex::Box& bx = mfi.validbox();
457
458 // Inputs
459 Set::Patch<const Set::Scalar> rho = rho_mf.array(mfi);
460 Set::Patch<const Set::Scalar> E = E_mf.array(mfi);
461 Set::Patch<const Set::Scalar> M = M_mf.array(mfi);
462
463 // Outputs
464 Set::Patch<Set::Scalar> rho_rhs = rho_rhs_mf.array(mfi);
465 Set::Patch<Set::Scalar> M_rhs = M_rhs_mf.array(mfi);
466 Set::Patch<Set::Scalar> E_rhs = E_rhs_mf.array(mfi);
467
468
469 // Set::Patch<Set::Scalar> rho_new = density_mf.Patch(lev,mfi);
470 // Set::Patch<Set::Scalar> E_new = energy_mf.Patch(lev,mfi);
471 // Set::Patch<Set::Scalar> M_new = momentum_mf.Patch(lev,mfi);
472
473 Set::Patch<const Set::Scalar> rho_solid = solid.density_mf.Patch(lev,mfi);
474 Set::Patch<const Set::Scalar> M_solid = solid.momentum_mf.Patch(lev,mfi);
475 Set::Patch<const Set::Scalar> E_solid = solid.energy_mf.Patch(lev,mfi);
476
477 Set::Patch<Set::Scalar> omega = vorticity_mf.Patch(lev,mfi);
478
479 Set::Patch<const Set::Scalar> eta = eta_old_mf.Patch(lev,mfi);
480 Set::Patch<const Set::Scalar> etadot = etadot_mf.Patch(lev,mfi);
481 Set::Patch<const Set::Scalar> velocity = velocity_mf.Patch(lev,mfi);
482
483 Set::Patch<const Set::Scalar> m0 = m0_mf.Patch(lev,mfi);
484 Set::Patch<const Set::Scalar> q = q_mf.Patch(lev,mfi);
485 Set::Patch<const Set::Scalar> _u0 = u0_mf.Patch(lev,mfi);
486
487 amrex::Array4<Set::Scalar> const& Source = (*Source_mf[lev]).array(mfi);
488
489 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
490 {
491 auto sten = Numeric::GetStencil(i, j, k, domain);
492
493 //Diffuse Sources
494 Set::Vector grad_eta = Numeric::Gradient(eta, i, j, k, 0, DX);
495 Set::Scalar grad_eta_mag = grad_eta.lpNorm<2>();
496 Set::Matrix hess_eta = Numeric::Hessian(eta, i, j, k, 0, DX);
497
498 Set::Vector u = Set::Vector(velocity(i, j, k, 0), velocity(i, j, k, 1));
499 Set::Vector u0 = Set::Vector(_u0(i, j, k, 0), _u0(i, j, k, 1));
500
501 Set::Matrix gradM = Numeric::Gradient(M, i, j, k, DX);
502 Set::Vector gradrho = Numeric::Gradient(rho,i,j,k,0,DX);
503 Set::Matrix hess_rho = Numeric::Hessian(rho,i,j,k,0,DX,sten);
504 Set::Matrix gradu = (gradM - u*gradrho.transpose()) / rho(i,j,k);
505
506 Set::Vector q0 = Set::Vector(q(i,j,k,0),q(i,j,k,1));
507
508 Set::Scalar mdot0 = -m0(i,j,k) * grad_eta_mag;
509 Set::Vector Pdot0 = Set::Vector::Zero();
510 Set::Scalar qdot0 = q0.dot(grad_eta);
511
512 // sten is necessary here because sometimes corner ghost
513 // cells don't get filled
514 Set::Matrix3 hess_M = Numeric::Hessian(M,i,j,k,DX);
516 for (int p = 0; p < 2; p++)
517 for (int q = 0; q < 2; q++)
518 for (int r = 0; r < 2; r++)
519 {
520 hess_u(r,p,q) =
521 (hess_M(r,p,q) - gradu(r,q)*gradrho(p) - gradu(r,p)*gradrho(q) - u(r)*hess_rho(p,q))
522 / rho(i,j,k);
523 }
524
525 Set::Vector Ldot0 = Set::Vector::Zero();
526 Set::Vector div_tau = Set::Vector::Zero();
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++)
531 {
532 Set::Scalar Mpqrs = 0.0;
533 if (p==r && q==s) Mpqrs += 0.5 * mu;
534
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);
537 }
538
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;// - Ldot0(0)*v(i,j,k,0) - Ldot0(1)*v(i,j,k,1);
543
544 // Lagrange terms to enforce no-penetration
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);
547
548 //Godunov flux
549 //states of total fields
550 const int X = 0, Y = 1;
551 Solver::Local::Riemann::State state_xlo(rho, M, E, i-1, j, k, X);
552 Solver::Local::Riemann::State state_x (rho, M, E, i , j, k, X);
553 Solver::Local::Riemann::State state_xhi(rho, M, E, i+1, j, k, X);
554
555 Solver::Local::Riemann::State state_ylo(rho, M, E, i, j-1, k, Y);
556 Solver::Local::Riemann::State state_y (rho, M, E, i, j , k, Y);
557 Solver::Local::Riemann::State state_yhi(rho, M, E, i, j+1, k, Y);
558
559 //states of solid fields
560 Solver::Local::Riemann::State state_xlo_solid(rho_solid, M_solid, E_solid, i-1, j, k, X);
561 Solver::Local::Riemann::State state_x_solid (rho_solid, M_solid, E_solid, i , j, k, X);
562 Solver::Local::Riemann::State state_xhi_solid(rho_solid, M_solid, E_solid, i+1, j, k, X);
563
564 Solver::Local::Riemann::State state_ylo_solid(rho_solid, M_solid, E_solid, i, j-1, k, Y);
565 Solver::Local::Riemann::State state_y_solid (rho_solid, M_solid, E_solid, i, j , k, Y);
566 Solver::Local::Riemann::State state_yhi_solid(rho_solid, M_solid, E_solid, i, j+1, k, Y);
567
568 Solver::Local::Riemann::State state_xlo_fluid = (state_xlo - (1.0 - eta(i-1,j,k))*state_xlo_solid) / (eta(i-1,j,k) + small);
569 Solver::Local::Riemann::State state_x_fluid = (state_x - (1.0 - eta(i,j,k) )*state_x_solid ) / (eta(i,j,k) + small);
570 Solver::Local::Riemann::State state_xhi_fluid = (state_xhi - (1.0 - eta(i+1,j,k))*state_xhi_solid) / (eta(i+1,j,k) + small);
571
572 Solver::Local::Riemann::State state_ylo_fluid = (state_ylo - (1.0 - eta(i,j-1,k))*state_ylo_solid) / (eta(i,j-1,k) + small);
573 Solver::Local::Riemann::State state_y_fluid = (state_y - (1.0 - eta(i,j,k) )*state_y_solid ) / (eta(i,j,k) + small);
574 Solver::Local::Riemann::State state_yhi_fluid = (state_yhi - (1.0 - eta(i,j+1,k))*state_yhi_solid) / (eta(i,j+1,k) + small);
575
576 Solver::Local::Riemann::Flux flux_xlo, flux_ylo, flux_xhi, flux_yhi;
577
578 try
579 {
580 //lo interface fluxes
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);
583
584 //hi interface fluxes
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);
587 }
588 catch(...)
589 {
590 Util::ParallelMessage(INFO,"lev=",lev);
591 Util::ParallelMessage(INFO,"i=",i,"j=",j);
593 }
594
595
596 Set::Scalar drhof_dt =
597 (flux_xlo.mass - flux_xhi.mass) / DX[0] +
598 (flux_ylo.mass - flux_yhi.mass) / DX[1] +
599 Source(i, j, k, 0);
600
601 rho_rhs(i,j,k) =
602 // rho_new(i, j, k) = rho(i, j, k) +
603 //(
604 drhof_dt +
605 // todo add drhos_dt term if want time-evolving rhos
606 etadot(i,j,k) * (rho(i,j,k) - rho_solid(i,j,k)) / (eta(i,j,k) + small)
607 // ) * dt;
608 ;
609
610
611
612 Set::Scalar dMxf_dt =
613 (flux_xlo.momentum_normal - flux_xhi.momentum_normal ) / DX[0] +
614 (flux_ylo.momentum_tangent - flux_yhi.momentum_tangent) / DX[1] +
615 div_tau(0) * eta(i,j,k) +
616 Source(i, j, k, 1);
617
618 M_rhs(i,j,k,0) =
619 //M_new(i, j, k, 0) = M(i, j, k, 0) +
620 // (
621 dMxf_dt +
622 // todo add dMs_dt term if want time-evolving Ms
623 etadot(i,j,k)*(M(i,j,k,0) - M_solid(i,j,k,0)) / (eta(i,j,k) + small)
624 // ) * dt;
625 ;
626
627 Set::Scalar dMyf_dt =
628 (flux_xlo.momentum_tangent - flux_xhi.momentum_tangent) / DX[0] +
629 (flux_ylo.momentum_normal - flux_yhi.momentum_normal ) / DX[1] +
630 div_tau(1) * eta(i,j,k) +
631 Source(i, j, k, 2);
632
633 M_rhs(i,j,k,1) =
634 //M_new(i, j, k, 1) = M(i, j, k, 1) +
635 //(
636 dMyf_dt +
637 // todo add dMs_dt term if want time-evolving Ms
638 etadot(i,j,k)*(M(i,j,k,1) - M_solid(i,j,k,1)) / (eta(i,j,k)+small)
639 // )*dt;
640 ;
641
642 Set::Scalar dEf_dt =
643 (flux_xlo.energy - flux_xhi.energy) / DX[0] +
644 (flux_ylo.energy - flux_yhi.energy) / DX[1] +
645 Source(i, j, k, 3);
646
647 E_rhs(i,j,k) =
648 // E_new(i, j, k) = E(i, j, k) +
649 // (
650 dEf_dt +
651 // todo add dEs_dt term if want time-evolving Es
652 etadot(i,j,k)*(E(i,j,k) - E_solid(i,j,k)) / (eta(i,j,k)+small)
653 // ) * dt;
654 ;
655
656#ifdef AMREX_DEBUG
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)))
661 {
662 Util::ParallelMessage(INFO,"rho_rhs=",rho_rhs(i,j,k));
663 Util::ParallelMessage(INFO,"Mx_rhs=",M_rhs(i,j,k,0));
664 Util::ParallelMessage(INFO,"Mx_rhs=",M_rhs(i,j,k,1));
665 Util::ParallelMessage(INFO,"E_rhs=",E_rhs(i,j,k));
666
667 Util::ParallelMessage(INFO,"lev=",lev);
668 Util::ParallelMessage(INFO,"i=",i," j=",j);
669 Util::ParallelMessage(INFO,"drhof_dt ",drhof_dt); // dies
670 Util::ParallelMessage(INFO,"flux_xlo.mass ",flux_xlo.mass);
671 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
672 Util::ParallelMessage(INFO,"flux_ylo.mass ",flux_ylo.mass);
673 Util::ParallelMessage(INFO,"flux_xhi.mass ",flux_yhi.mass);
674 Util::ParallelMessage(INFO,"eta ",eta(i,j,k));
675 Util::ParallelMessage(INFO,"etadot ",etadot(i,j,k));
676 Util::ParallelMessage(INFO,"Source ",Source(i,j,k,0));
677 Util::ParallelMessage(INFO,"state_x ",state_x); // <<<<
678 Util::ParallelMessage(INFO,"state_y ",state_y);
679 Util::ParallelMessage(INFO,"state_x_solid ",state_x_solid); // <<<<
680 Util::ParallelMessage(INFO,"state_y_solid ",state_y_solid);
681 Util::ParallelMessage(INFO,"state_xhi ",state_xhi); // <<<<
682 Util::ParallelMessage(INFO,"state_yhi ",state_yhi);
683 Util::ParallelMessage(INFO,"state_xhi_solid ",state_xhi_solid);
684 Util::ParallelMessage(INFO,"state_yhi_solids ",state_yhi_solid);
685 Util::ParallelMessage(INFO,"state_xlo ",state_xlo);
686 Util::ParallelMessage(INFO,"state_ylo ",state_ylo);
687 Util::ParallelMessage(INFO,"state_xlo_solid ",state_xlo_solid);
688 Util::ParallelMessage(INFO,"state_ylo_solid ",state_ylo_solid);
689
690 Util::ParallelMessage(INFO,"Mx_solid ",M_solid(i,j,k,0));
691 Util::ParallelMessage(INFO,"My_solid ",M_solid(i,j,k,1));
692 Util::ParallelMessage(INFO,"small ",small);
693 Util::ParallelMessage(INFO,"Mx ",M(i,j,k,0));
694 Util::ParallelMessage(INFO,"My ",M(i,j,k,1));
695 Util::ParallelMessage(INFO,"dMx/dt ",dMxf_dt);
696 Util::ParallelMessage(INFO,"dMy/dt ",dMyf_dt);
697
698
701 Util::Message(INFO,DX[0]);
704 Util::Message(INFO,DX[1]);
705 Util::Message(INFO,div_tau);
706 Util::Message(INFO,Source(i, j, k, 2));
707
708 Util::Message(INFO,hess_eta);
709 Util::Message(INFO,velocity(i,j,k,0));
710 Util::Message(INFO,velocity(i,j,k,1));
711
713 }
714#endif
715
716
717
718 // todo - may need to move this for higher order schemes...
719 omega(i, j, k) = eta(i, j, k) * (gradu(1,0) - gradu(0,1));
720 });
721 }
722}
723
724void Hydro::Regrid(int lev, Set::Scalar /* time */)
725{
726 BL_PROFILE("Integrator::Hydro::Regrid");
727 Source_mf[lev]->setVal(0.0);
728 if (lev < finest_level) return;
729
730 Util::Message(INFO, "Regridding on level", lev);
731}//end regrid
732
733//void Hydro::TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar time, int ngrow)
734void Hydro::TagCellsForRefinement(int lev, amrex::TagBoxArray& a_tags, Set::Scalar, int)
735{
736 BL_PROFILE("Integrator::Flame::TagCellsForRefinement");
737
738 const Set::Scalar* DX = geom[lev].CellSize();
739 Set::Scalar dr = sqrt(AMREX_D_TERM(DX[0] * DX[0], +DX[1] * DX[1], +DX[2] * DX[2]));
740
741 // Eta criterion for refinement
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);
746
747 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
748 Set::Vector grad_eta = Numeric::Gradient(eta, i, j, k, 0, DX);
749 if (grad_eta.lpNorm<2>() * dr * 2 > eta_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
750 });
751 }
752
753 // Vorticity criterion for refinement
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);
758
759 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
760 auto sten = Numeric::GetStencil(i, j, k, bx);
761 Set::Vector grad_omega = Numeric::Gradient(omega, i, j, k, 0, DX, sten);
762 if (grad_omega.lpNorm<2>() * dr * 2 > omega_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
763 });
764 }
765
766 // Gradu criterion for refinement
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);
771
772 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
773 auto sten = Numeric::GetStencil(i, j, k, bx);
774 Set::Matrix grad_u = Numeric::Gradient(v, i, j, k, DX, sten);
775 if (grad_u.lpNorm<2>() * dr * 2 > gradu_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
776 });
777 }
778
779 // Pressure criterion for refinement
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);
784
785 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
786 auto sten = Numeric::GetStencil(i, j, k, bx);
787 Set::Vector grad_p = Numeric::Gradient(p, i, j, k, 0, DX, sten);
788 if (grad_p.lpNorm<2>() * dr * 2 > p_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
789 });
790 }
791
792 // Density criterion for refinement
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);
797
798 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
799 auto sten = Numeric::GetStencil(i, j, k, bx);
800 Set::Vector grad_rho = Numeric::Gradient(rho, i, j, k, 0, DX, sten);
801 if (grad_rho.lpNorm<2>() * dr * 2 > rho_refinement_criterion) tags(i, j, k) = amrex::TagBox::SET;
802 });
803 }
804
805}
806
807}
808
809
810#endif
#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
Definition BMP.H:22
Initialize Laminates in a matrix.
Definition Laminate.H:16
Definition PNG.H:26
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
static Matrix3 Zero()
Definition Matrix3.H:42
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