Alamo
Flame.cpp
Go to the documentation of this file.
1#include "Flame.H"
2#include "IO/ParmParse.H"
3#include "BC/Constant.H"
4#include "Numeric/Stencil.H"
5#include "IC/Laminate.H"
6#include "IC/Constant.H"
7#include "IC/PSRead.H"
8#include "Numeric/Function.H"
9#include "IC/Expression.H"
10#include "IC/BMP.H"
11#include "IC/PNG.H"
12#include "Base/Mechanics.H"
13#include "Util/Util.H"
17#include <cmath>
18
19namespace Integrator
20{
21
24
26{
27 pp_queryclass(*this);
28}
29
30
31void
33{
34 pp.forbid("pressure.P","use chamber.pressure instead");
35
36 pp.forbid("geometry.x_len","This is specified by geometry.prob_lo/hi");
37 pp.forbid("geometry.y_len","This is specified by geometry.prob_lo/hi");
38 pp.forbid("amr.ghost_cells", "This should not be adjustable ");
39
40 pp.forbid("pf.gamma","use propellant.powerlaw.gamma");
41
42 pp.forbid("pressure.r_ap", "use propellant.powerlaw.r_ap");
43 pp.forbid("pressure.r_htpb", "use propellant.powerlaw.r_htpb");
44 pp.forbid("pressure.r_comb", "use propellant.powerlaw.r_comb");
45 pp.forbid("pressure.n_ap", "use propellant.powerlaw.n_ap");
46 pp.forbid("pressure.n_htpb", "use propellant.powerlaw.n_htpb");
47 pp.forbid("pressure.n_comb", "use propellant.powerlaw.n_comb");
48
49 pp.forbid("thermal.bound", "use thermal.Tref");
50 pp.forbid("thermal.T_fluid", "use thermal.Tfluid (or nothing)");
51 pp.forbid("thermal.m_ap", "use propellant.fullfeedback.m_ap");
52 pp.forbid("thermal.m_htpb", "use propellant.fullfeedback.m_htpb");
53 pp.forbid("thermal.E_ap", "use propellant.fullfeedback.E_ap");
54 pp.forbid("thermal.E_htpb", "use propellant.fullfeedback.E_htpb");
55 pp.forbid("thermal.modeling_ap", "Old debug variable. Should equal 1 ");
56 pp.forbid("thermal.modeling_htpb", "Old debug variable. Should equal 1");
57
58 pp.forbid("pressure.a1", "use propellant.fullfeedback.a1 instead");
59 pp.forbid("pressure.a2", "use propellant.fullfeedback.a2 instead");
60 pp.forbid("pressure.a3", "use propellant.fullfeedback.a3 instead");
61 pp.forbid("pressure.b1", "use propellant.fullfeedback.b1 instead");
62 pp.forbid("pressure.b2", "use propellant.fullfeedback.b2 instead");
63 pp.forbid("pressure.b3", "use propellant.fullfeedback.b3 instead");
64 pp.forbid("pressure.c1", "use propellant.fullfeedback.c1 instead");
65 pp.forbid("pressure.mob_ap", "no longer used");
66 pp.forbid("pressure.dependency", "use propellant.fullfeedback.pressure_dependency");
67 pp.forbid("pressure.h1", "use propellant.homogenize.h1 instead");
68 pp.forbid("pressure.h2", "use propellant.homogenize.h2 instead");
69 pp.forbid("thermal.mlocal_ap", "use propellant.homogenize.mlocal_ap");
70 pp.forbid("thermal.mlocal_comb", "use propellant.homogenize.mlocal_comb");
71 pp.forbid("thermal.mlocal_htpb", "this actually did **nothing** - it was overridden by a hard code using massfraction.");
72
73 pp.forbid("thermal.disperssion1", "use propellant.homogenize.dispersion1");
74 pp.forbid("thermal.disperssion2", "use propellant.homogenize.dispersion2");
75 pp.forbid("thermal.disperssion3", "use propellant.homogenize.dispersion3");
76
77 pp.forbid("thermal.rho_ap", "use propellant.fullfeedback/homogenize.rho_ap ");
78 pp.forbid("thermal.rho_htpb","use propellant.fullfeedback/homogenize.rho_htpb ");
79 pp.forbid("thermal.k_ap", "use propellant.fullfeedback/homogenize.k_ap ");
80 pp.forbid("thermal.k_htpb", "use propellant.fullfeedback/homogenize.k_htpb ");
81 pp.forbid("thermal.cp_ap", "use propellant.fullfeedback/homogenize.cp_ap ");
82 pp.forbid("thermal.cp_htpb","use propellant.fullfeedback/homogenize.cp_htpb ");
83}
84
85
86// [parser]
87void
89{
90 BL_PROFILE("Integrator::Flame::Flame()");
91
92 Forbids(pp);
93
94 // Whether to include extra fields (such as mdot, etc) in the plot output
95 pp.query_default("plot_field",value.plot_field,true);
96
97 //
98 // PHASE FIELD VARIABLES
99 //
100
101 // Burn width thickness
102 pp.query_default("pf.eps", value.pf.eps, "1.0_m", Unit::Length());
103 // Interface energy param
104 pp.query_default("pf.kappa", value.pf.kappa, "0.0_J/m^2", Unit::Energy() / Unit::Area());
105 // Chemical potential multiplier
106 pp.query_default("pf.lambda", value.pf.lambda, "0.0_J/m^2", Unit::Energy()/Unit::Area());
107 // Unburned rest energy
108 pp.query_default("pf.w1", value.pf.w1, "0.0",Unit::Less());
109 // Barrier energy
110 pp.query_default("pf.w12", value.pf.w12, "0.0", Unit::Less());
111 // Burned rest energy
112 pp.query_default("pf.w0", value.pf.w0, "0.0",Unit::Less());
113
114 // Boundary conditions for phase field order params
115 pp.select<BC::Constant>("pf.eta.bc", value.bc_eta, 1 );
116 value.RegisterNewFab(value.eta_mf, value.bc_eta, 1, 2, "eta", true);
117 value.RegisterNewFab(value.eta_old_mf, value.bc_eta, 1, 2, "eta_old", 0);
118
119 // Inital value of eta that doesn't evolve and is used during refiment to set the updated values of eta with voids in the domain.
120 // Used to fix a bug where duirn refinement, a void won't be updated correctly and would be a square, not a circle
121 value.RegisterNewFab(value.eta_0_mf, value.bc_eta, 1, 2, "eta_0", 0);
122
123 // value.RegisterNewFab(value.eta_mf_frozen, value.bc_eta, 1, 2, "eta_frozen", value.plot_field);
124
125 // phase field initial condition
127
128
129 // Select reduced order model to capture heat feedback
133 ("propellant",value.propellant);
134
135
136 // Whether to use the Thermal Transport Model
137 pp_query_default("thermal.on", value.thermal.on, false);
138
139 // Reference temperature
140 // Used to set all other reference temperatures by default.
141 pp_query_default("thermal.Tref", value.thermal.Tref, "300.0_K",Unit::Temperature());
142
143 if (value.thermal.on) {
144
145 // Used to change heat flux units
146 pp_query_default("thermal.hc", value.thermal.hc, "1.0", Unit::Power()/Unit::Area());
147
148 // Effective fluid temperature, temp of the eta = 0 (fluid) region
149 pp_query_default("thermal.Tfluid", value.thermal.Tfluid, value.thermal.Tref);
150
151 // Cutoff value for regression, if T < Tcutoff eta won't evolve/regress
152 pp.query_default("thermal.Tcutoff", value.thermal.Tcutoff, "0.0", Unit::Temperature());
153
154 // Switch time of the improved regridding where eta and the temperature field are both used. It is recommended to make this time ~10x the timestep
155 // Before this the refinement is based on the gradient of eta which helps the laser IC start correctly. A regrid is forced when this time is reached.
156 pp.query_default("thermal.end_initial_refine_time", value.thermal.end_initial_refine_time, "0.0", Unit::Time());
157
158 //Temperature boundary condition
159 pp.select_default<BC::Constant>("thermal.temp.bc", value.bc_temp, 1, Unit::Temperature());
160
161 value.RegisterNewFab(value.temp_mf, value.bc_temp, 1, 3, "temp", true);
162 value.RegisterNewFab(value.temp_old_mf, value.bc_temp, 1, 3, "temp_old", false);
163 value.RegisterNewFab(value.temps_mf, value.bc_temp, 1, 0, "temps", false);
164
165 value.RegisterNewFab(value.mdot_mf, value.bc_temp, 1, 0, "mdot", value.plot_field);
166 value.RegisterNewFab(value.alpha_mf, value.bc_temp, 1, 0, "alpha", value.plot_field);
167 value.RegisterNewFab(value.heatflux_mf, value.bc_temp, 1, 0, "heatflux", value.plot_field);
168 value.RegisterNewFab(value.laser_mf, value.bc_temp, 1, 0, "laser", value.plot_field);
169
170 value.RegisterIntegratedVariable(&value.chamber.volume, "volume");
171 value.RegisterIntegratedVariable(&value.chamber.area, "area");
172 value.RegisterIntegratedVariable(&value.chamber.massflux, "mass_flux");
173
174 value.RegisterNewFab(value.thermal.has_exceeded_Tcutoff, value.bc_temp, 1, 2, "exceeded_Tcutoff", 0); // Used to determine where regression has started
175
176 // laser initial condition
179 ("laser.ic",value.ic_laser, value.geom, Unit::Power()/Unit::Area());
180
181 // thermal initial condition
184 IC::BMP,
185 IC::PNG >
186 ("temp.ic",value.thermal.ic_temp,value.geom, Unit::Temperature());
187 }
188
189
190 // Constant pressure value
191 pp_query_default("chamber.pressure", value.chamber.pressure, "1.0_MPa", Unit::Pressure());
192
193 // Whether to compute the pressure evolution
194 pp_query_default("variable_pressure", value.variable_pressure, false);
195
196 // Refinement criterion for eta field
197 pp_query_default( "amr.refinement_criterion", value.m_refinement_criterion, "0.001",
198 Unit::Less());
199
200 // Refinement criterion for temperature field
201 pp.query_default( "amr.refinement_criterion_temp", value.t_refinement_criterion, "0.001_K",
203
204 // Eta value to restrict the refinament for the temperature field
205 pp.query_default( "amr.refinament_restriction", value.t_refinement_restriction, "0.1",
206 Unit::Less());
207
208 // Refinement criterion for phi field [infinity]
209 pp_query_default("amr.phi_refinement_criterion", value.phi_refinement_criterion, 1.0e100);
210
211 // Minimum allowable threshold for $\eta$
212 pp_query_default("small", value.small, 1.0e-8);
213
214 // Initial condition for $\phi$ field.
216 ("phi.ic",value.ic_phi,value.geom);
217
218 value.RegisterNodalFab(value.phi_mf, 1, 2, "phi", true);
219
220 // Whether to use Neo-hookean Elastic model
221 pp_query_default("elastic.on", value.elastic.on, 0);
222
223 // Body force
224 pp_query_default("elastic.traction", value.elastic.traction, 0.0);
225
226 // Phi refinement criteria
227 pp_query_default("elastic.phirefinement", value.elastic.phirefinement, 1);
228
229 pp.queryclass<Base::Mechanics<model_type>>("elastic",value);
230
231 if (value.m_type != Type::Disable)
232 {
233 // Reference temperature for thermal expansion
234 // (temperature at which the material is strain-free)
235 pp_query_default("Telastic", value.elastic.Telastic, value.thermal.Tref);
236 // elastic model of AP
238 // elastic model of HTPB
240
241 // Use our current eta field as the psi field for the solver
242 value.psi_on = false;
243 value.solver.setPsi(value.eta_mf);
244 }
245
246 bool allow_unused;
247 // Set this to true to allow unused inputs without error.
248 // (Not recommended.)
249 pp.query_default("allow_unused",allow_unused,false);
250 if (!allow_unused && pp.AnyUnusedInputs())
251 {
252 Util::Warning(INFO,"The following inputs were specified but not used:");
253 pp.AllUnusedInputs();
254 Util::Exception(INFO,"Aborting. Specify 'allow_unused=True` to ignore this error.");
255 }
256}
257
258void Flame::Initialize(int lev)
259{
260 BL_PROFILE("Integrator::Flame::Initialize");
262
263 ic_eta->Initialize(lev, eta_mf);
265 ic_phi->Initialize(lev, phi_mf);
266 //ic_phicell->Initialize(lev, phicell_mf);
267
268 if (elastic.on) {
269 rhs_mf[lev]->setVal(Set::Vector::Zero());
270 }
271 if (thermal.on) {
272 if (thermal.ic_temp)
273 {
274 thermal.ic_temp->Initialize(lev,temp_mf);
275 thermal.ic_temp->Initialize(lev,temp_old_mf);
276 thermal.ic_temp->Initialize(lev,temps_mf);
277 }
278 else
279 {
280 temp_mf[lev]->setVal(thermal.Tref);
281 temp_old_mf[lev]->setVal(thermal.Tref);
282 temps_mf[lev]->setVal(thermal.Tref);
283 }
284 alpha_mf[lev]->setVal(0.0);
285 mdot_mf[lev]->setVal(0.0);
286 heatflux_mf[lev]->setVal(0.0);
288 }
289 if (variable_pressure) chamber.pressure = 1.0;
290}
291
292void Flame::UpdateModel(int /*a_step*/, Set::Scalar /*a_time*/)
293{
295
296 for (int lev = 0; lev <= finest_level; ++lev)
297 {
298 amrex::Box domain = this->geom[lev].Domain();
299 domain.convert(amrex::IntVect::TheNodeVector());
300 const Set::Scalar* DX = geom[lev].CellSize();
301
302 phi_mf[lev]->FillBoundary();
303 eta_mf[lev]->FillBoundary();
304 temp_mf[lev]->FillBoundary();
305
306 for (MFIter mfi(*model_mf[lev], false); mfi.isValid(); ++mfi)
307 {
308 amrex::Box smallbox = mfi.nodaltilebox();
309 amrex::Box bx = mfi.grownnodaltilebox() & domain;
310 Set::Patch<model_type> model = model_mf.Patch(lev,mfi);
313 Set::Patch<Set::Vector> rhs = rhs_mf.Patch(lev,mfi);
314
315 if (elastic.on)
316 {
318 amrex::ParallelFor(smallbox, [=] AMREX_GPU_DEVICE(int i, int j, int k)
319
320 {
321 Set::Vector grad_eta = Numeric::CellGradientOnNode(eta, i, j, k, 0, DX);
322
323 rhs(i, j, k) = (elastic.traction) * grad_eta;
324
325 });
326 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
327 {
328 Set::Scalar phi_avg = phi(i, j, k, 0);
329 Set::Scalar temp_avg = Numeric::Interpolate::CellToNodeAverage(temp, i, j, k, 0);
330 model_type model_ap = elastic.model_ap;
331 model_ap.F0 -= Set::Matrix::Identity();
332 model_ap.F0 *= (temp_avg - elastic.Telastic);
333 model_ap.F0 += Set::Matrix::Identity();
334 model_type model_htpb = elastic.model_htpb;
335 model_htpb.F0 -= Set::Matrix::Identity();
336 model_htpb.F0 *= (temp_avg - elastic.Telastic);
337 model_htpb.F0 += Set::Matrix::Identity();
338
339 model(i, j, k) = (model_ap * phi_avg + model_htpb * (1. - phi_avg));
340 });
341 }
342 else
343 {
344 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
345 {
346 Set::Scalar phi_avg = Numeric::Interpolate::CellToNodeAverage(phi, i, j, k, 0);
347 //phi_avg = phi(i,j,k,0);
348 model_type model_ap = elastic.model_ap;
349 model_ap.F0 *= Set::Matrix::Zero();
350 model_type model_htpb = elastic.model_htpb;
351 model_htpb.F0 *= Set::Matrix::Zero();
352 model(i, j, k) = (model_ap * phi_avg + model_htpb * (1. - phi_avg));
353 });
354 }
355 }
356 Util::RealFillBoundary(*model_mf[lev], geom[lev]);
357
358 }
359}
360
361void Flame::TimeStepBegin(Set::Scalar a_time, int a_iter)
362{
363 BL_PROFILE("Integrator::Flame::TimeStepBegin");
365 if (thermal.on) {
366 for (int lev = 0; lev <= finest_level; ++lev)
367 ic_laser->Initialize(lev, laser_mf, a_time);
368 }
369
370 if (a_time > thermal.end_initial_refine_time)
371 {
372 if (!end_initial_refine) {
373 for (int lev = 0; lev <= finest_level; ++lev)
374 Flame::Regrid(lev, a_time);
376 }
377
378 prev_finest_ba = grids[finest_level];
379 prev_finest_level = finest_level;
380 }
381}
382
383void Flame::TimeStepComplete(Set::Scalar /*a_time*/, int /*a_iter*/)
384{
385 BL_PROFILE("Integrator::Flame::TimeStepComplete");
386 if (variable_pressure) {
387 //Set::Scalar x_len = geom[0].ProbDomain().length(0);
388 //Set::Scalar y_len = geom[0].ProbDomain().length(1);
389 // Set::Scalar domain_area = x_len * y_len;
390 Util::Message(INFO, "Mass = ", chamber.massflux);
391 Util::Message(INFO, "Pressure = ", chamber.pressure);
392 }
393}
394
396{
397 BL_PROFILE("Integrador::Flame::Advance");
399 const Set::Scalar* DX = geom[lev].CellSize();
400
401 std::swap(eta_old_mf[lev], eta_mf[lev]);
402
403 //
404 // Chamber pressure update
405 //
406 if (variable_pressure) {
407 chamber.pressure = exp(0.00075 * chamber.massflux);
408 if (chamber.pressure > 10.0) {
409 chamber.pressure = 10.0;
410 }
411 else if (chamber.pressure <= 0.99) {
412 chamber.pressure = 0.99;
413 }
414 elastic.traction = chamber.pressure;
415 }
416
417
418 //
419 // Multi-well chemical potential
420 //
422 0.0,
423 -5.0 * pf.w1 + 16.0 * pf.w12 - 11.0 * pf.w0,
424 14.0 * pf.w1 - 32.0 * pf.w12 + 18.0 * pf.w0,
425 -8.0 * pf.w1 + 16.0 * pf.w12 - 8.0 * pf.w0);
427
429
430 for (amrex::MFIter mfi(*eta_mf[lev], true); mfi.isValid(); ++mfi)
431 {
432 const amrex::Box& bx = mfi.tilebox();
433 // Phase fields
434 Set::Patch<Set::Scalar> etanew = eta_mf.Patch(lev,mfi);
437 // Heat transfer fields
439 Set::Patch<Set::Scalar> alpha = alpha_mf.Patch(lev,mfi);
440 Set::Patch<Set::Scalar> laser = laser_mf.Patch(lev,mfi);
441 // Diagnostic fields
442 Set::Patch<Set::Scalar> mdot = mdot_mf.Patch(lev,mfi);
443 Set::Patch<Set::Scalar> heatflux = heatflux_mf.Patch(lev,mfi);
444
445 Set::Patch<Set::Scalar> exceeded_Tcutoff = thermal.has_exceeded_Tcutoff.Patch(lev, mfi);
446 Set::Scalar Tcutoff = thermal.Tcutoff;
447
448
449 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
450 {
451 //
452 // CALCULATE PHI-AVERAGED QUANTITIES
453 //
454 Set::Scalar phi_avg = Numeric::Interpolate::NodeToCellAverage(phi, i, j, k, 0);
455 Set::Scalar T = thermal.on ? temp(i,j,k) : NAN;
456
457 Set::Scalar K = propellant.get_K(phi_avg);
458
460
461 Set::Scalar cp = propellant.get_cp(phi_avg);
462
463 //
464 // CALCULATE MOBILITY
465 //
466 Set::Scalar L = propellant.get_L( phi_avg, T);
467
468 //
469 // EVOLVE PHASE FIELD (ETA)
470 //
471
472 Set::Scalar eta_lap = Numeric::Laplacian(eta, i, j, k, 0, DX);
473 Set::Scalar df_deta = ((pf.lambda / pf.eps) * dw(eta(i, j, k)) - pf.eps * pf.kappa * eta_lap);
474
475 if (df_deta < 0) {
476 // Prevent eta from increasing/healing. A bug was found where if the diffuse thickness was too large compared to a void
477 // (region of eta = 0), eta would heal/increase in a non-physcial way, this statement stops that behavior
478 df_deta = 0.0;
479 }
480 if (thermal.on && T < thermal.Tcutoff) {
481 // If the temperature is lower then the cutoff temperature don't evolve the eta field
482 df_deta = 0.0;
483 }
484 etanew(i, j, k) = eta(i, j, k) - L * dt * df_deta;
485
486 if (etanew(i, j, k) <= small) etanew(i, j, k) = small;
487
488 if (thermal.on)
489 {
490 //
491 // Calculate thermal diffisivity and store for later gradient
492 //
493
494 alpha(i, j, k) = K / rho / cp;
495
496 //
497 // CALCULATE MASS FLUX BASED ON EVOLVING ETA
498 //
499
500 mdot(i, j, k) = rho * fabs(eta(i, j, k) - etanew(i, j, k)) / dt;
501
502 //
503 // CALCULATE HEAT FLUX BASED ON THE CALCULATED MASS FLUX
504 //
505
506 Set::Scalar q0 = propellant.get_qdot(mdot(i,j,k), phi_avg);
507 heatflux(i,j,k) = ( thermal.hc*q0 + laser(i,j,k) ) / K;
508
509 if (temp(i,j,k) > Tcutoff)
510 {
511 exceeded_Tcutoff(i,j,k) = 1;
512 }
513
514 }
515
516 });
517
518 } // MFi For loop
519
520
521 //
522 // THERMAL TRANSPORT
523 //
524 if (thermal.on)
525 {
526 std::swap(temp_old_mf[lev], temp_mf[lev]);
527
528 for (amrex::MFIter mfi(*eta_mf[lev], true); mfi.isValid(); ++mfi)
529 {
530 const amrex::Box& bx = mfi.tilebox();
531
532 Set::Patch<Set::Scalar> tempnew = temp_mf.Patch(lev,mfi);
535
536 Set::Patch<Set::Scalar> temps = temps_mf.Patch(lev,mfi);
537
538
539 // Phase field
540 Set::Patch<Set::Scalar> etanew = (*eta_mf[lev]).array(mfi);
541 Set::Patch<const Set::Scalar> eta = (*eta_old_mf[lev]).array(mfi);
542 // Diagnostic fields
544
545 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
546 {
547 auto sten = Numeric::GetStencil(i, j, k, bx);
548 Set::Vector grad_eta = Numeric::Gradient(eta, i, j, k, 0, DX);
549 Set::Vector grad_temp = Numeric::Gradient(temp, i, j, k, 0, DX);
550 Set::Scalar lap_temp = Numeric::Laplacian(temp, i, j, k, 0, DX);
551 Set::Scalar grad_eta_mag = grad_eta.lpNorm<2>();
552 Set::Vector grad_alpha = Numeric::Gradient(alpha, i, j, k, 0, DX, sten);
553 Set::Scalar dTdt = 0.0;
554 dTdt += grad_eta.dot(grad_temp * alpha(i, j, k));
555 dTdt += grad_alpha.dot(eta(i, j, k) * grad_temp);
556 dTdt += eta(i, j, k) * alpha(i, j, k) * lap_temp;
557 dTdt += alpha(i, j, k) * heatflux(i, j, k) * grad_eta_mag;
558
559 Set::Scalar Tsolid = dTdt + temps(i, j, k) * (etanew(i, j, k) - eta(i, j, k)) / dt;
560 temps(i, j, k) = temps(i, j, k) + dt * Tsolid;
561 tempnew(i, j, k) = etanew(i, j, k) * temps(i, j, k) + (1.0 - etanew(i, j, k)) * thermal.Tfluid;
562 });
563 }
564 }
565
566} //Function
567
568
569void Flame::TagCellsForRefinement(int lev, amrex::TagBoxArray& a_tags, Set::Scalar time, int ngrow)
570{
571 BL_PROFILE("Integrator::Flame::TagCellsForRefinement");
573
574 const Set::Scalar* DX = geom[lev].CellSize();
575 Set::Scalar dr = sqrt(AMREX_D_TERM(DX[0] * DX[0], +DX[1] * DX[1], +DX[2] * DX[2]));
576
577 // Eta criterion for refinement
578 for (amrex::MFIter mfi(*eta_mf[lev], true); mfi.isValid(); ++mfi)
579 {
580 const amrex::Box& bx = mfi.tilebox();
581 amrex::Array4<char> const& tags = a_tags.array(mfi);
584
585 if (thermal.on) {
586 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
587 {
588 Set::Vector gradeta = Numeric::Gradient(eta, i, j, k, 0, DX);
589 if (gradeta.lpNorm<2>() * dr * 2 > m_refinement_criterion && eta(i, j, k) >= t_refinement_restriction && temp(i,j,k) > thermal.Tcutoff*0.9)
590 tags(i, j, k) = amrex::TagBox::SET;
591 });
592
593 } else {
594 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
595 {
596 Set::Vector gradeta = Numeric::Gradient(eta, i, j, k, 0, DX);
597 if (gradeta.lpNorm<2>() * dr * 2 > m_refinement_criterion && eta(i, j, k) >= t_refinement_restriction)
598 tags(i, j, k) = amrex::TagBox::SET;
599 });
600 }
601 }
602
603 // Phi criterion for refinement
604 if (elastic.phirefinement) {
605 for (amrex::MFIter mfi(*eta_mf[lev], true); mfi.isValid(); ++mfi)
606 {
607 const amrex::Box& bx = mfi.tilebox();
608 amrex::Array4<char> const& tags = a_tags.array(mfi);
610
611 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
612 {
613 Set::Vector gradphi = Numeric::Gradient(phi, i, j, k, 0, DX);
614 if (gradphi.lpNorm<2>() * dr >= phi_refinement_criterion)
615 tags(i, j, k) = amrex::TagBox::SET;
616 });
617 }
618 }
619
620
621 // Thermal criterion for refinement
622 if (thermal.on) {
623 for (amrex::MFIter mfi(*temp_mf[lev], true); mfi.isValid(); ++mfi)
624 {
625 const amrex::Box& bx = mfi.tilebox();
626 amrex::Array4<char> const& tags = a_tags.array(mfi);
629 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
630 {
631 Set::Vector tempgrad = Numeric::Gradient(temp, i, j, k, 0, DX);
632 if (tempgrad.lpNorm<2>() * dr > t_refinement_criterion && eta(i, j, k) >= t_refinement_restriction)
633 tags(i, j, k) = amrex::TagBox::SET;
634 });
635 }
636 }
637
638 // Refine at start
639 for (amrex::MFIter mfi(*eta_mf[lev], true); mfi.isValid(); ++mfi)
640 {
641 const amrex::Box& bx = mfi.tilebox();
642 amrex::Array4<char> const& tags = a_tags.array(mfi);
644
645 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
646 {
647 Set::Vector gradeta = Numeric::Gradient(eta, i, j, k, 0, DX);
648 if (gradeta.lpNorm<2>() * dr * 2 > m_refinement_criterion && time < thermal.end_initial_refine_time)
649 tags(i, j, k) = amrex::TagBox::SET;
650 });
651 }
652}
653
654void Flame::Regrid(int lev, Set::Scalar time)
655{
656 BL_PROFILE("Integrator::Flame::Regrid");
657
658 ic_phi->Initialize(lev, phi_mf, time);
659 ic_eta->Initialize(lev, eta_0_mf, time);
660
661 if (thermal.on) {
662 /*
663 This regrid function works by making a using the "has_exceeded_Tcutoff" field. If the temperature in a cell is greater than Tcutoff,
664 eta will change and when regruding won't use the initial eta field. If T < T_cutoff, when regriding happens it applies the inital
665 eta field condition. This gives at leat a 4x speed improvement in 2D when doing regression with voids. This is because orgioanlly
666 there was a bug where when regridding, the orgional eta field wouldn't be applied, so there would be "sqaures" of voids instead of
667 circles/spheres when using .xyzr files as the inital condition.
668 */
669 for (amrex::MFIter mfi(*eta_mf[lev], true); mfi.isValid(); ++mfi)
670 {
671 const amrex::Box &bx = mfi.tilebox();
672 Set::Patch<Set::Scalar> eta = eta_mf.Patch(lev, mfi);
675 Set::Patch<Set::Scalar> exceeded_Tcutoff = thermal.has_exceeded_Tcutoff.Patch(lev, mfi);
676
677 Set::Scalar Tcutoff = thermal.Tcutoff;
678
679 amrex::BoxList boxes_to_update;
680 if (lev == finest_level && prev_finest_level == finest_level)
681 boxes_to_update = amrex::complementIn(bx, prev_finest_ba).boxList();
682 else
683 boxes_to_update.push_back(bx);
684
685 for (const amrex::Box &box : boxes_to_update)
686 amrex::ParallelFor(box, [=] AMREX_GPU_DEVICE(int i, int j, int k)
687 {
688
689 if (!exceeded_Tcutoff(i,j,k) && temp(i,j,k) < Tcutoff)
690 {
691 eta(i, j, k) = eta_0(i, j, k);
692 }
693 });
694 }
695
696 if (lev == finest_level)
697 {
698 prev_finest_ba = grids[finest_level];
699 prev_finest_level = finest_level;
700 }
701 }
702}
703
704void Flame::Integrate(int amrlev, Set::Scalar time, int /*step*/,
705 const amrex::MFIter& mfi, const amrex::Box& box)
706{
707 BL_PROFILE("Flame::Integrate");
708
710
711 const Set::Scalar* DX = geom[amrlev].CellSize();
712 Set::Scalar dv = AMREX_D_TERM(DX[0], *DX[1], *DX[2]);
715 if (variable_pressure) {
716 amrex::ParallelFor(box, [=] AMREX_GPU_DEVICE(int i, int j, int k)
717 {
718 chamber.volume += eta(i, j, k, 0) * dv;
719 Set::Vector grad = Numeric::Gradient(eta, i, j, k, 0, DX);
720 Set::Scalar normgrad = grad.lpNorm<2>();
721 Set::Scalar da = normgrad * dv;
722 chamber.area += da;
723
724 Set::Vector mgrad = Numeric::Gradient(mdot, i, j, k, 0, DX);
725 Set::Scalar mnormgrad = mgrad.lpNorm<2>();
726 Set::Scalar dm = mnormgrad * dv;
727 chamber.massflux += dm;
728
729 });
730 }
731 else {
732 amrex::ParallelFor(box, [=] AMREX_GPU_DEVICE(int i, int j, int k)
733 {
734 chamber.volume += eta(i, j, k, 0) * dv;
735 Set::Vector grad = Numeric::Gradient(eta, i, j, k, 0, DX);
736 Set::Scalar normgrad = grad.lpNorm<2>();
737 Set::Scalar da = normgrad * dv;
738 chamber.area += da;
739 });
740 }
741 // time dependent pressure data from experimenta -> p = 0.0954521220950523 * exp(15.289993148880678 * t)
742}
743} // namespace Integrator
744
#define pp_query_default(...)
Definition ParmParse.H:102
#define pp_queryclass(...)
Definition ParmParse.H:109
#define INFO
Definition Util.H:24
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
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(std::string name, PTRTYPE *&ic_eta, Args &&... args)
Definition ParmParse.H:1056
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
virtual void TimeStepBegin(Set::Scalar a_time, int a_step) override
Definition Mechanics.H:170
void Integrate(int amrlev, Set::Scalar, int, const amrex::MFIter &mfi, const amrex::Box &a_box) override
Definition Mechanics.H:391
void Initialize(int lev) override
Use the #ic object to initialize::Temp.
Definition Mechanics.H:150
void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
Definition Mechanics.H:251
Solver::Nonlocal::Newton< MODEL > solver
Definition Mechanics.H:506
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar, int) override
Definition Mechanics.H:450
Set::Field< Model::Solid::Finite::NeoHookeanPredeformed > model_mf
Definition Mechanics.H:474
Set::Scalar Tcutoff
Definition Flame.H:108
Model::Propellant::Propellant< Model::Propellant::PowerLaw, Model::Propellant::FullFeedback, Model::Propellant::Homogenize > propellant
Definition Flame.H:138
Set::Scalar volume
Definition Flame.H:126
Set::Scalar Telastic
Definition Flame.H:117
IC::IC< Set::Scalar > * ic_eta
Definition Flame.H:85
Set::Field< Set::Scalar > laser_mf
Definition Flame.H:73
Set::Field< Set::Scalar > alpha_mf
Definition Flame.H:70
amrex::BoxArray prev_finest_ba
Definition Flame.H:89
int prev_finest_level
Definition Flame.H:90
struct Integrator::Flame::@3 elastic
BC::BC< Set::Scalar > * bc_temp
Definition Flame.H:75
void UpdateModel(int a_step, Set::Scalar a_time) override
Definition Flame.cpp:292
Set::Scalar lambda
Definition Flame.H:96
Set::Scalar w1
Definition Flame.H:98
Set::Field< Set::Scalar > eta_0_mf
Definition Flame.H:65
Set::Field< Set::Scalar > phi_mf
Definition Flame.H:68
Set::Field< Set::Scalar > temp_old_mf
Definition Flame.H:60
Set::Scalar small
Definition Flame.H:84
IC::IC< Set::Scalar > * ic_temp
Definition Flame.H:111
int variable_pressure
Definition Flame.H:87
Set::Scalar area
Definition Flame.H:127
void Regrid(int lev, Set::Scalar time) override
Definition Flame.cpp:654
Set::Scalar pressure
Definition Flame.H:130
model_type model_ap
Definition Flame.H:118
Set::Scalar w0
Definition Flame.H:98
void TimeStepComplete(Set::Scalar a_time, int a_iter) override
Definition Flame.cpp:383
static void Forbids(IO::ParmParse &pp)
Definition Flame.cpp:32
struct Integrator::Flame::@2 thermal
Set::Field< Set::Scalar > eta_mf
Definition Flame.H:63
Set::Scalar phi_refinement_criterion
Definition Flame.H:80
void TimeStepBegin(Set::Scalar a_time, int a_iter) override
Definition Flame.cpp:361
static void Parse(Flame &value, IO::ParmParse &pp)
Definition Flame.cpp:88
void Initialize(int lev) override
Definition Flame.cpp:258
Set::Scalar massflux
Definition Flame.H:128
struct Integrator::Flame::@1 pf
Set::Scalar kappa
Definition Flame.H:97
Set::Scalar end_initial_refine_time
Definition Flame.H:109
Set::Scalar traction
Definition Flame.H:119
Set::Scalar Tref
Definition Flame.H:106
Set::Field< Set::Scalar > temps_mf
Definition Flame.H:61
void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
Definition Flame.cpp:395
Set::Scalar t_refinement_restriction
Definition Flame.H:83
struct Integrator::Flame::@4 chamber
Set::Scalar t_refinement_criterion
Definition Flame.H:82
Set::Field< Set::Scalar > eta_old_mf
Definition Flame.H:64
Set::Scalar m_refinement_criterion
Definition Flame.H:81
Set::Field< Set::Scalar > has_exceeded_Tcutoff
Definition Flame.H:110
Set::Scalar hc
Definition Flame.H:105
void Integrate(int amrlev, Set::Scalar time, int step, const amrex::MFIter &mfi, const amrex::Box &box) override
Definition Flame.cpp:704
model_type model_htpb
Definition Flame.H:118
Set::Scalar Tfluid
Definition Flame.H:107
BC::BC< Set::Scalar > * bc_eta
Definition Flame.H:76
int end_initial_refine
Definition Flame.H:91
IC::IC< Set::Scalar > * ic_laser
Definition Flame.H:78
void TagCellsForRefinement(int lev, amrex::TagBoxArray &tags, amrex::Real, int) override
Definition Flame.cpp:569
Set::Field< Set::Scalar > temp_mf
Definition Flame.H:59
Set::Scalar w12
Definition Flame.H:98
Set::Field< Set::Scalar > heatflux_mf
Definition Flame.H:71
bool plot_field
Definition Flame.H:86
Set::Field< Set::Scalar > mdot_mf
Definition Flame.H:66
IC::IC< Set::Scalar > * ic_phi
Definition Flame.H:77
Set::Scalar eps
Definition Flame.H:95
amrex::Real timestep
Timestep for the base level of refinement.
Definition Integrator.H:393
void RegisterNodalFab(Set::Field< Set::Scalar > &new_fab, int ncomp, int nghost, std::string name, bool writeout, bool evolving=true, std::vector< std::string > suffix={})
Add a new node-based scalar field.
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.
std::vector< amrex::Box > box
Definition Integrator.H:465
amrex::Vector< amrex::Real > dt
Timesteps for each level of refinement.
Definition Integrator.H:394
void RegisterIntegratedVariable(Set::Scalar *integrated_variable, std::string name, bool extensive=true)
Register a variable to be integrated over the spatial domain using the Integrate function.
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Set::Scalar get_cp(const Set::Scalar phi)
Definition Propellant.H:78
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Set::Scalar get_qdot(const Set::Scalar mdot, const Set::Scalar phi)
Definition Propellant.H:91
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE void set_pressure(Set::Scalar P)
Definition Propellant.H:40
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Set::Scalar get_L(const Set::Scalar phi, const Set::Scalar T)
Definition Propellant.H:104
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Set::Scalar get_rho(const Set::Scalar phi)
Definition Propellant.H:65
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Set::Scalar get_K(const Set::Scalar phi)
Definition Propellant.H:52
amrex::Array4< Set::Scalar > Patch(const int lev, const amrex::MFIter &mfi) const &
Definition Set.H:263
amrex::Array4< T > Patch(int lev, amrex::MFIter &mfi) const &
Definition Set.H:76
void setPsi(Set::Field< Set::Scalar > &a_psi)
Definition Newton.H:46
Collection of numerical integrator objects.
Definition AllenCahn.H:42
AMREX_FORCE_INLINE Set::Scalar Laplacian(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:555
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::Vector CellGradientOnNode(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])
Definition Stencil.H:699
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
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
AMREX_FORCE_INLINE void RealFillBoundary(amrex::FabArray< amrex::BaseFab< T > > &a_mf, const amrex::Geometry &, const int nghost=2)
Definition Util.H:339
static AMREX_FORCE_INLINE T NodeToCellAverage(const amrex::Array4< const T > &f, const int &i, const int &j, const int &k, const int &m)
Definition Stencil.H:1428
static AMREX_FORCE_INLINE T CellToNodeAverage(const amrex::Array4< const T > &f, const int &i, const int &j, const int &k, const int &m, std::array< StencilType, AMREX_SPACEDIM > stencil=DefaultType)
Definition Stencil.H:1390
static Unit Area()
Definition Unit.H:209
static Unit Energy()
Definition Unit.H:218
static Unit Time()
Definition Unit.H:199
static Unit Temperature()
Definition Unit.H:201
static Unit Pressure()
Definition Unit.H:217
static Unit Power()
Definition Unit.H:219
static Unit Length()
Definition Unit.H:198
static Unit Less()
Definition Unit.H:197