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 "Base/Mechanics.H"
12 
13 #include <cmath>
14 
15 namespace Integrator
16 {
17 
19 
21 {
22  pp_queryclass(*this);
23 }
24 
25 // [parser]
26 void
28 {
29  BL_PROFILE("Integrator::Flame::Flame()");
30  {
31  //pp_query("fields_verbose", value.plot_field);
32  pp_query_default("timestep", value.base_time, 1.0e-4);
33  // These are the phase field method parameters
34  // that you use to inform the phase field method.
35  pp_query_default("pf.eps", value.pf.eps, 0.0); // Burn width thickness
36  pp_query_default("pf.kappa", value.pf.kappa, 0.0); // Interface energy param
37  pp_query_default("pf.gamma", value.pf.gamma, 1.0); // Scaling factor for mobility
38  pp_query_default("pf.lambda", value.pf.lambda, 0.0); // Chemical potential multiplier
39  pp_query_default("pf.w1", value.pf.w1, 0.0); // Unburned rest energy
40  pp_query_default("pf.w12", value.pf.w12, 0.0); // Barrier energy
41  pp_query_default("pf.w0", value.pf.w0, 0.0); // Burned rest energy
42  pp_query_default("amr.ghost_cells", value.ghost_count, 2); // number of ghost cells in all fields
43  pp_query_default("geometry.x_len", value.x_len, 0.001); // Domain x length
44  pp_query_default("geometry.y_len", value.y_len, 0.001); // Domain y length
45 
46  value.bc_eta = new BC::Constant(1);
47  pp_queryclass("pf.eta.bc", *static_cast<BC::Constant*>(value.bc_eta)); // See :ref:`BC::Constant`
48  value.RegisterNewFab(value.eta_mf, value.bc_eta, 1, value.ghost_count, "eta", true);
49  value.RegisterNewFab(value.eta_old_mf, value.bc_eta, 1, value.ghost_count, "eta_old", false);
50 
51  std::string eta_bc_str = "constant";
52  pp_query_validate("pf.eta.ic.type", eta_bc_str, {"constant", "expression"}); // Eta boundary condition [constant, expression]
53  if (eta_bc_str == "constant") value.ic_eta = new IC::Constant(value.geom, pp, "pf.eta.ic.constant");
54  else if (eta_bc_str == "expression") value.ic_eta = new IC::Expression(value.geom, pp, "pf.eta.ic.expression");
55 
56  std::string eta_ic_type = "constant";
57  pp_query_validate("eta.ic.type", eta_ic_type, {"laminate", "constant", "expression"}); // Eta initial condition [constant, laminate, expression, bmp]
58  if (eta_ic_type == "laminate") value.ic_eta = new IC::Laminate(value.geom, pp, "eta.ic.laminate");
59  else if (eta_ic_type == "constant") value.ic_eta = new IC::Constant(value.geom, pp, "eta.ic.constant");
60  else if (eta_ic_type == "expression") value.ic_eta = new IC::Expression(value.geom, pp, "eta.ic.expression");
61  else if (eta_ic_type == "bmp") value.ic_eta = new IC::BMP(value.geom, pp, "eta.ic.bmp");
62  else if (eta_ic_type == "png") value.ic_eta = new IC::PNG(value.geom, pp, "eta.ic.png");
63  else Util::Abort(INFO, "Invalid eta IC type", eta_ic_type);
64  }
65 
66  {
67  //IO::ParmParse pp("thermal");
68  pp_query_default("thermal.on", value.thermal.on, false); // Whether to use the Thermal Transport Model
69  pp_query_default("elastic.on", value.elastic.on, 0); // Whether to use Neo-hookean Elastic model
70  pp_query_default("thermal.bound", value.thermal.bound, 0.0); // System Initial Temperature
71  pp_query_default("elastic.traction", value.elastic.traction, 0.0); // Body force
72  pp_query_default("elastic.phirefinement", value.elastic.phirefinement, 1); // Phi refinement criteria
73 
74 
75 
76  if (value.thermal.on) {
77  pp_query_required("thermal.rho_ap", value.thermal.rho_ap); // AP Density
78  pp_query_required("thermal.rho_htpb", value.thermal.rho_htpb); // HTPB Density
79  pp_query_required("thermal.k_ap", value.thermal.k_ap); // AP Thermal Conductivity
80  pp_query_required("thermal.k_htpb", value.thermal.k_htpb); // HTPB Thermal Conductivity
81  pp_query_required("thermal.cp_ap", value.thermal.cp_ap); // AP Specific Heat
82  pp_query_required("thermal.cp_htpb", value.thermal.cp_htpb); //HTPB Specific Heat
83 
84  pp_query_default("thermal.q0", value.thermal.q0, 0.0); // Baseline heat flux
85 
86  pp_query_required("thermal.m_ap", value.thermal.m_ap); // AP Pre-exponential factor for Arrhenius Law
87  pp_query_required("thermal.m_htpb", value.thermal.m_htpb); // HTPB Pre-exponential factor for Arrhenius Law
88  pp_query_required("thermal.E_ap", value.thermal.E_ap); // AP Activation Energy for Arrhenius Law
89  pp_query_required("thermal.E_htpb", value.thermal.E_htpb); // HTPB Activation Energy for Arrhenius Law
90 
91  pp_query_default("thermal.hc", value.thermal.hc, 1.0); // Used to change heat flux units
92  pp_query_default("thermal.massfraction", value.thermal.massfraction, 0.8); // Systen AP mass fraction
93  pp_query_default("thermal.mlocal_ap", value.thermal.mlocal_ap, 0.0); // AP mass flux reference value
94  pp_query_default("thermal.mlocal_htpb", value.thermal.mlocal_htpb, 0.0); // HTPB mass flux reference value
95  pp_query_default("thermal.mlocal_comb", value.thermal.mlocal_comb, 0.0); // AP/HTPB mass flux reference value
96 
97  pp_query_default("thermal.T_fluid", value.thermal.T_fluid, 300.0); // Temperature of the Standin Fluid
98 
99  pp_query("thermal.disperssion1", value.thermal.disperssion1); // K; dispersion variables are use to set the outter field properties for the void grain case.
100  pp_query("thermal.disperssion2", value.thermal.disperssion2); // rho; dispersion variables are use to set the outter field properties for the void grain case.
101  pp_query("thermal.disperssion3", value.thermal.disperssion3); // cp; dispersion variables are use to set the outter field properties for the void grain case.
102 
103  pp_query_default("thermal.modeling_ap", value.thermal.modeling_ap, 1.0); // Scaling factor for AP thermal conductivity (default = 1.0)
104  pp_query_default("thermal.modeling_htpb", value.thermal.modeling_htpb, 1.0); // Scaling factor for HTPB thermal conductivity (default = 1.0)
105 
106 
107  value.bc_temp = new BC::Constant(1);
108  pp_queryclass("thermal.temp.bc", *static_cast<BC::Constant*>(value.bc_temp));
109  value.RegisterNewFab(value.temp_mf, value.bc_temp, 1, value.ghost_count + 1, "temp", true);
110  value.RegisterNewFab(value.temp_old_mf, value.bc_temp, 1, value.ghost_count + 1, "temp_old", false);
111  value.RegisterNewFab(value.temps_mf, value.bc_temp, 1, value.ghost_count + 1, "temps", false);
112  value.RegisterNewFab(value.temps_old_mf, value.bc_temp, 1, value.ghost_count + 1, "temps_old", false);
113 
114  value.RegisterNewFab(value.mdot_mf, value.bc_eta, 1, value.ghost_count + 1, "mdot", value.plot_field);
115  value.RegisterNewFab(value.mob_mf, value.bc_eta, 1, value.ghost_count + 1, "mob", value.plot_field);
116  value.RegisterNewFab(value.alpha_mf, value.bc_temp, 1, value.ghost_count + 1, "alpha", value.plot_field);
117  value.RegisterNewFab(value.heatflux_mf, value.bc_temp, 1, value.ghost_count + 1, "heatflux", value.plot_field);
118  value.RegisterNewFab(value.laser_mf, value.bc_temp, 1, value.ghost_count + 1, "laser", value.plot_field);
119 
120  value.RegisterIntegratedVariable(&value.volume, "total_area");
121  value.RegisterIntegratedVariable(&value.area, "Interface_area");
122  value.RegisterIntegratedVariable(&value.chamber_area, "chamber_area", false);
123  value.RegisterIntegratedVariable(&value.massflux, "mass_flux");
124  value.RegisterIntegratedVariable(&value.chamber_pressure, "Pressure", false);
125 
126  std::string laser_ic_type = "constant";
127  pp_query_validate("laser.ic.type", laser_ic_type, {"expression", "constant"}); // heat laser initial condition type [constant, expression]
128  if (laser_ic_type == "expression") value.ic_laser = new IC::Expression(value.geom, pp, "laser.ic.expression");
129  else if (laser_ic_type == "constant") value.ic_laser = new IC::Constant(value.geom, pp, "laser.ic.constant");
130  else Util::Abort(INFO, "Invalid eta IC type", laser_ic_type);
131 
132  std::string temp_ic_type;
133  pp_query_validate("temp.ic.type", temp_ic_type,{"default","constant","expression","bmp","png"}); // Temperature initial condition
134  if (temp_ic_type == "constant") value.thermal.ic_temp = new IC::Constant(value.geom, pp, "temp.ic.constant");
135  else if (temp_ic_type == "expression") value.thermal.ic_temp = new IC::Expression(value.geom, pp, "temp.ic.expression");
136  else if (temp_ic_type == "bmp") value.thermal.ic_temp = new IC::BMP(value.geom, pp, "temp.ic.bmp");
137  else if (temp_ic_type == "png") value.thermal.ic_temp = new IC::PNG(value.geom, pp, "temp.ic.png");
138  else if (temp_ic_type == "default") value.thermal.ic_temp = nullptr;
139  }
140  }
141 
142  {
143  pp_query_default("pressure.P", value.pressure.P, 1.0); // Constant pressure value
144  if (value.thermal.on)
145  {
146  pp_query_required("pressure.a1", value.pressure.arrhenius.a1); // Surgate heat flux model paramater - AP
147  pp_query_required("pressure.a2", value.pressure.arrhenius.a2); // Surgate heat flux model paramater - HTPB
148  pp_query_required("pressure.a3", value.pressure.arrhenius.a3); // Surgate heat flux model paramater - Total
149  pp_query_required("pressure.b1", value.pressure.arrhenius.b1); // Surgate heat flux model paramater - AP
150  pp_query_required("pressure.b2", value.pressure.arrhenius.b2); // Surgate heat flux model paramater - HTPB
151  pp_query_required("pressure.b3", value.pressure.arrhenius.b3); // Surgate heat flux model paramater - Total
152  pp_query_required("pressure.c1", value.pressure.arrhenius.c1); // Surgate heat flux model paramater - Total
153  pp_query_default("pressure.mob_ap", value.pressure.arrhenius.mob_ap, 0); // Whether to include pressure to the arrhenius law
154  pp_query_default("pressure.dependency", value.pressure.arrhenius.dependency, 1); // Whether to use pressure to determined the reference Zeta
155  pp_query_default("pressure.h1", value.pressure.arrhenius.h1, 1.81); // Surgate heat flux model paramater - Homogenized
156  pp_query_default("pressure.h2", value.pressure.arrhenius.h2, 1.34); // Surgate heat flux model paramater - Homogenized
157 
158  }
159  else
160  {
161  pp_query_required("pressure.r_ap", value.pressure.power.r_ap); // AP power pressure law parameter (r*P^n)
162  pp_query_required("pressure.r_htpb", value.pressure.power.r_htpb); // HTPB power pressure law parameter (r*P^n)
163  pp_query_required("pressure.r_comb", value.pressure.power.r_comb); // AP/HTPB power pressure law parameter (r*P^n)
164  pp_query_required("pressure.n_ap", value.pressure.power.n_ap); // AP power pressure law parameter (r*P^n)
165  pp_query_required("pressure.n_htpb", value.pressure.power.n_htpb); // HTPB power pressure law parameter (r*P^n)
166  pp_query_required("pressure.n_comb", value.pressure.power.n_comb); // AP/HTPB power pressure law parameter (r*P^n)
167  }
168  pp_query_default("variable_pressure", value.variable_pressure, 0); // Whether to compute the pressure evolution
169  pp_query_default("homogeneousSystem", value.homogeneousSystem, 0); // Whether to initialize Phi with homogenized properties
170  }
171 
172  pp_query_default("amr.refinement_criterion", value.m_refinement_criterion, 0.001);// Refinement criterion for eta field
173  pp_query_default("amr.refinement_criterion_temp", value.t_refinement_criterion, 0.001);// Refinement criterion for temperature field
174  pp_query_default("amr.refinament_restriction", value.t_refinement_restriction, 0.1);// Eta value to restrict the refinament for the temperature field
175  pp_query_default("amr.phi_refinement_criterion", value.phi_refinement_criterion, 1.0e100);// Refinement criterion for phi field [infinity]
176  pp_query_default("small", value.small, 1.0e-8); // Lowest value of Eta.
177 
178  {
179  // The material field is referred to as :math:`\phi(\mathbf{x})` and is
180  // specified using these parameters.
181  //IO::ParmParse pp("phi.ic");
182  std::string phi_ic_type = "packedspheres";
183  pp_query_validate("phi.ic.type", phi_ic_type, {"psread", "laminate", "expression", "constant", "bmp", "png"}); // IC type (psread, laminate, constant)
184  if (phi_ic_type == "psread") {
185  value.ic_phi = new IC::PSRead(value.geom, pp, "phi.ic.psread");
186  //value.ic_phicell = new IC::PSRead(value.geom, pp, "phi.ic.psread");
187  pp_query_default("phi.ic.psread.eps", value.zeta, 1.0e-5); // AP/HTPB interface length
188  pp_query_default("phi.zeta_0", value.zeta_0, 1.0e-5); // Reference interface length for heat integration
189  }
190  else if (phi_ic_type == "laminate") {
191  value.ic_phi = new IC::Laminate(value.geom, pp, "phi.ic.laminate");
192  //value.ic_phicell = new IC::Laminate(value.geom, pp, "phi.ic.laminate");
193  pp_query_default("phi.ic.laminate.eps", value.zeta, 1.0e-5); // AP/HTPB interface length
194  pp_query_default("phi.zeta_0", value.zeta_0, 1.0e-5); // Reference interface length for heat integration
195  }
196  else if (phi_ic_type == "expression") {
197  value.ic_phi = new IC::Expression(value.geom, pp, "phi.ic.expression");
198  //value.ic_phicell = new IC::Expression(value.geom, pp, "phi.ic.expression");
199  pp_query_default("phi.zeta_0", value.zeta_0, 1.0e-5); // Reference interface length for heat integration
200  pp_query_default("phi.zeta", value.zeta, 1.0e-5); // AP/HTPB interface length
201  }
202  else if (phi_ic_type == "constant") {
203  value.ic_phi = new IC::Constant(value.geom, pp, "phi.ic.constant");
204  //value.ic_phicell = new IC::Constant(value.geom, pp, "phi.ic.constant");
205  pp_query_default("phi.zeta_0", value.zeta_0, 1.0e-5); // Reference interface length for heat integration
206  pp_query_default("phi.zeta", value.zeta, 1.0e-5); // AP/HTPB interface length
207  }
208  else if (phi_ic_type == "bmp") {
209  value.ic_phi = new IC::BMP(value.geom, pp, "phi.ic.bmp");
210  //value.ic_phicell = new IC::BMP(value.geom, pp, "phi.ic.bmp");
211  pp_query_default("phi.zeta_0", value.zeta_0, 1.0e-5); // Reference interface length for heat integration
212  pp_query_default("phi.zeta", value.zeta, 1.0e-5); // AP/HTPB interface length
213  }
214  else if (phi_ic_type == "png") {
215  value.ic_phi = new IC::PNG(value.geom, pp, "phi.ic.png");
216  pp_query_default("phi.zeta_0", value.zeta_0, 1.0e-5); // Reference interface length for heat integration
217  pp_query_default("phi.zeta", value.zeta, 1.0e-5); // AP/HTPB interface length
218  }
219  else Util::Abort(INFO, "Invalid IC type ", phi_ic_type);
220  //value.RegisterNewFab(value.phi_mf, 1, "phi_cell", true);
221  value.RegisterNodalFab(value.phi_mf, 1, value.ghost_count + 1, "phi", true);
222  //value.RegisterNewFab(value.phicell_mf, value.bc_eta, 1, value.ghost_count + 1, "phi", true);
223  }
224 
225  pp_queryclass("elastic",static_cast<Base::Mechanics<model_type>&>(value));
226 
227  if (value.m_type != Type::Disable)
228  {
229  value.elastic.Tref = value.thermal.bound;
230  pp_query_default("Tref", value.elastic.Tref, 300.0); // Initial temperature for thermal expansion computation
231  pp_queryclass("model_ap", value.elastic.model_ap);
232  pp_queryclass("model_htpb", value.elastic.model_htpb);
233 
234  value.bc_psi = new BC::Nothing();
235  value.RegisterNewFab(value.psi_mf, value.bc_psi, 1, value.ghost_count, "psi", value.plot_psi);
236  value.psi_on = true;
237  }
238 }
239 
240 void Flame::Initialize(int lev)
241 {
242  BL_PROFILE("Integrator::Flame::Initialize");
245 
246  ic_eta->Initialize(lev, eta_mf);
248  ic_phi->Initialize(lev, phi_mf);
249  //ic_phicell->Initialize(lev, phicell_mf);
250 
251  if (elastic.on) {
252  psi_mf[lev]->setVal(1.0);
253  rhs_mf[lev]->setVal(Set::Vector::Zero());
254  }
255  if (thermal.on) {
256  if (thermal.ic_temp)
257  {
258  thermal.ic_temp->Initialize(lev,temp_mf);
259  thermal.ic_temp->Initialize(lev,temp_old_mf);
260  thermal.ic_temp->Initialize(lev,temps_mf);
261  thermal.ic_temp->Initialize(lev,temps_old_mf);
262  }
263  else
264  {
265  temp_mf[lev]->setVal(thermal.bound);
266  temp_old_mf[lev]->setVal(thermal.bound);
267  temps_mf[lev]->setVal(thermal.bound);
268  temps_old_mf[lev]->setVal(thermal.bound);
269  }
270  alpha_mf[lev]->setVal(0.0);
271  mob_mf[lev]->setVal(0.0);
272  mdot_mf[lev]->setVal(0.0);
273  heatflux_mf[lev]->setVal(0.0);
274  // pressure_mf[lev]->setVal(1.0); // [error]
275  thermal.w1 = 0.2 * pressure.P + 0.9;
276  thermal.T_fluid = thermal.bound;
278  thermal.mlocal_htpb = 685000.0 - 850e3 * thermal.massfraction;
279  }
280  if (variable_pressure) pressure.P = 1.0;
281 }
282 
283 void Flame::UpdateModel(int /*a_step*/, Set::Scalar /*a_time*/)
284 {
286 
287  for (int lev = 0; lev <= finest_level; ++lev)
288  {
289  amrex::Box domain = this->geom[lev].Domain();
290  domain.convert(amrex::IntVect::TheNodeVector());
291  const Set::Scalar* DX = geom[lev].CellSize();
292 
293  //psi_mf[lev]->setVal(1.0);
294  phi_mf[lev]->FillBoundary();
295  //phicell_mf[lev]->FillBoundary();
296  eta_mf[lev]->FillBoundary();
297  temp_mf[lev]->FillBoundary();
298 
299  for (MFIter mfi(*model_mf[lev], false); mfi.isValid(); ++mfi)
300  {
301  amrex::Box bx = mfi.grownnodaltilebox() & domain;
302  //amrex::Box bx = mfi.nodaltilebox();
303  //bx.grow(1);
304  amrex::Array4<model_type> const& model = model_mf[lev]->array(mfi);
305  amrex::Array4<const Set::Scalar> const& phi = phi_mf[lev]->array(mfi);
306  amrex::Array4<const Set::Scalar> const& eta = eta_mf[lev]->array(mfi);
307  amrex::Array4<Set::Vector> const& rhs = rhs_mf[lev]->array(mfi);
308  // amrex::Array4<const Set::Scalar> const& Pressure = pressure_mf[lev]->array(mfi); // [error]
309 
310  if (elastic.on)
311  {
312  amrex::Array4<const Set::Scalar> const& temp = temp_mf[lev]->array(mfi);
313  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
314  {
315  //auto sten = Numeric::GetStencil(i, j, k, bx);
316  Set::Scalar phi_avg = phi(i, j, k, 0);
317  Set::Scalar temp_avg = Numeric::Interpolate::CellToNodeAverage(temp, i, j, k, 0);
318  Set::Vector grad_eta = Numeric::CellGradientOnNode(eta, i, j, k, 0, DX);
319  rhs(i, j, k) = elastic.traction * grad_eta;
320 
321  model_type model_ap = elastic.model_ap;
322  model_ap.F0 -= Set::Matrix::Identity();
323  model_ap.F0 *= (temp_avg - elastic.Tref);
324  model_ap.F0 += Set::Matrix::Identity();
325  model_type model_htpb = elastic.model_htpb;
326  model_htpb.F0 -= Set::Matrix::Identity();
327  model_htpb.F0 *= (temp_avg - elastic.Tref);
328  model_htpb.F0 += Set::Matrix::Identity();
329 
330  model(i, j, k) = model_ap * phi_avg + model_htpb * (1. - phi_avg);
331  });
332  }
333  else
334  {
335  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
336  {
337  Set::Scalar phi_avg = Numeric::Interpolate::CellToNodeAverage(phi, i, j, k, 0);
338  //phi_avg = phi(i,j,k,0);
339  model_type model_ap = elastic.model_ap;
340  model_ap.F0 *= Set::Matrix::Zero();
341  model_type model_htpb = elastic.model_htpb;
342  model_htpb.F0 *= Set::Matrix::Zero();
343  model(i, j, k) = model_ap * phi_avg + model_htpb * (1. - phi_avg);
344  });
345  }
346  }
347  Util::RealFillBoundary(*model_mf[lev], geom[lev]);
348 
349  psi_mf[lev]->setVal(1.0);
350  amrex::MultiFab::Copy(*psi_mf[lev], *eta_mf[lev], 0, 0, 1, psi_mf[lev]->nGrow());
351  }
352 }
353 
354 void Flame::TimeStepBegin(Set::Scalar a_time, int a_iter)
355 {
356  BL_PROFILE("Integrator::Flame::TimeStepBegin");
358  if (thermal.on) {
359  for (int lev = 0; lev <= finest_level; ++lev)
360  ic_laser->Initialize(lev, laser_mf, a_time);
361  }
362 }
363 
364 void Flame::TimeStepComplete(Set::Scalar /*a_time*/, int /*a_iter*/)
365 {
366  BL_PROFILE("Integrator::Flame::TimeStepComplete");
367  if (variable_pressure) {
368  Set::Scalar domain_area = x_len * y_len;
370  chamber_area = domain_area - volume;
371  Util::Message(INFO, "Mass = ", massflux);
372  Util::Message(INFO, "Pressure = ", pressure.P);
373  }
374 }
375 
376 void Flame::Advance(int lev, Set::Scalar time, Set::Scalar dt)
377 {
378  BL_PROFILE("Integrador::Flame::Advance");
380  const Set::Scalar* DX = geom[lev].CellSize();
381 
382  if (true) //lev == finest_level) //(true)
383  {
384  if (thermal.on)
385  {
386  std::swap(eta_old_mf[lev], eta_mf[lev]);
387  std::swap(temp_old_mf[lev], temp_mf[lev]);
388  std::swap(temps_old_mf[lev], temps_mf[lev]);
389 
391  0.0,
392  -5.0 * pf.w1 + 16.0 * pf.w12 - 11.0 * pf.w0,
393  14.0 * pf.w1 - 32.0 * pf.w12 + 18.0 * pf.w0,
394  -8.0 * pf.w1 + 16.0 * pf.w12 - 8.0 * pf.w0);
396 
397  for (amrex::MFIter mfi(*eta_mf[lev], true); mfi.isValid(); ++mfi)
398  {
399  const amrex::Box& bx = mfi.tilebox();
400  // Phase fields
401  amrex::Array4<Set::Scalar> const& etanew = (*eta_mf[lev]).array(mfi);
402  amrex::Array4<const Set::Scalar> const& eta = (*eta_old_mf[lev]).array(mfi);
403  amrex::Array4<const Set::Scalar> const& phi = (*phi_mf[lev]).array(mfi);
404  //amrex::Array4<const Set::Scalar> const& phicell = (*phicell_mf[lev]).array(mfi);
405  // Heat transfer fields
406  amrex::Array4<Set::Scalar> const& tempnew = (*temp_mf[lev]).array(mfi);
407  amrex::Array4<const Set::Scalar> const& temp = (*temp_old_mf[lev]).array(mfi);
408  amrex::Array4<Set::Scalar> const& alpha = (*alpha_mf[lev]).array(mfi);
409  amrex::Array4<Set::Scalar> const& tempsnew = (*temps_mf[lev]).array(mfi);
410  amrex::Array4<Set::Scalar> const& temps = (*temps_old_mf[lev]).array(mfi);
411  amrex::Array4<Set::Scalar> const& laser = (*laser_mf[lev]).array(mfi);
412  // Diagnostic fields
413  amrex::Array4<Set::Scalar> const& mob = (*mob_mf[lev]).array(mfi);
414  amrex::Array4<Set::Scalar> const& mdot = (*mdot_mf[lev]).array(mfi);
415  amrex::Array4<Set::Scalar> const& heatflux = (*heatflux_mf[lev]).array(mfi);
416 
417  if (variable_pressure) {
418  pressure.P = exp(0.00075 * massflux);
419  if (pressure.P > 10.0) {
420  pressure.P = 10.0;
421  }
422  else if (pressure.P <= 0.99) {
423  pressure.P = 0.99;
424  }
425  elastic.traction = pressure.P;
426 
427  }
428 
429  Set::Scalar zeta_2 = 0.000045 - pressure.P * 6.42e-6;
430  Set::Scalar zeta_1;
431  if (pressure.arrhenius.dependency == 1) zeta_1 = zeta_2;
432  else zeta_1 = zeta_0;
433 
434  Set::Scalar k1 = pressure.arrhenius.a1 * pressure.P + pressure.arrhenius.b1 - zeta_1 / zeta;
435  Set::Scalar k2 = pressure.arrhenius.a2 * pressure.P + pressure.arrhenius.b2 - zeta_1 / zeta;
436  Set::Scalar k3 = 4.0 * log((pressure.arrhenius.c1 * pressure.P * pressure.P + pressure.arrhenius.a3 * pressure.P + pressure.arrhenius.b3) - k1 / 2.0 - k2 / 2.0);
437  Set::Scalar k4 = pressure.arrhenius.h1 * pressure.P + pressure.arrhenius.h2;
438 
439  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
440  {
441  Set::Scalar phi_avg = Numeric::Interpolate::NodeToCellAverage(phi, i, j, k, 0);
442  Set::Scalar eta_lap = Numeric::Laplacian(eta, i, j, k, 0, DX);
443  Set::Scalar K; // Calculate effective thermal conductivity
444  Set::Scalar rho; // No special interface mixure rule is needed here.
445  Set::Scalar cp;
446  if (homogeneousSystem) {
447  K = (thermal.modeling_ap * thermal.k_ap * thermal.massfraction + thermal.modeling_htpb * thermal.k_htpb * (1.0 - thermal.massfraction)) * phi_avg + thermal.disperssion1 * (1. - phi_avg); // Calculate effective thermal conductivity
448  rho = (thermal.rho_ap * thermal.massfraction + thermal.rho_htpb * (1.0 - thermal.massfraction)) * phi_avg + thermal.disperssion2 * (1. - phi_avg); // No special interface mixure rule is needed here.
449  cp = (thermal.cp_ap * thermal.massfraction + thermal.cp_htpb * (1.0 - thermal.massfraction)) * phi_avg + thermal.disperssion3 * (1. - phi_avg);
450  }
451  else {
452  K = thermal.modeling_ap * thermal.k_ap * phi_avg + thermal.modeling_htpb * thermal.k_htpb * (1.0 - phi_avg); // Calculate effective thermal conductivity
453  rho = thermal.rho_ap * phi_avg + thermal.rho_htpb * (1.0 - phi_avg); // No special interface mixure rule is needed here.
454  cp = thermal.cp_ap * phi_avg + thermal.cp_htpb * (1.0 - phi_avg);
455  }
456  Set::Scalar df_deta = ((pf.lambda / pf.eps) * dw(eta(i, j, k)) - pf.eps * pf.kappa * eta_lap);
457  etanew(i, j, k) = eta(i, j, k) - mob(i, j, k) * dt * df_deta;
458  if (etanew(i, j, k) <= small) etanew(i, j, k) = small;
459 
460  alpha(i, j, k) = K / rho / cp; // Calculate thermal diffusivity and store in fiel
461  mdot(i, j, k) = rho * fabs(eta(i, j, k) - etanew(i, j, k)) / dt; // deta/dt
462 
463  if (isnan(etanew(i, j, k)) || isnan(alpha(i, j, k)) || isnan(mdot(i, j, k))) {
464  Util::Message(INFO, etanew(i, j, k), "etanew contains nan (i=", i, " j= ", j, ")");
465  Util::Message(INFO, mdot(i, j, k), "mdot contains nan (i=", i, " j= ", j, ")");
466  Util::Message(INFO, alpha(i, j, k), "alpha contains nan (i=", i, " j= ", j, ")");
467  Util::Abort(INFO);
468  }
469  });
470 
471  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
472  {
473  Set::Scalar phi_avg = Numeric::Interpolate::NodeToCellAverage(phi, i, j, k, 0);
474  //phi_avg = phicell(i,j,k);
475  Set::Scalar K;
476  Set::Scalar qflux;
477  Set::Scalar mlocal;
478  Set::Scalar mdota;
479  Set::Scalar mbase;
480 
481  if (homogeneousSystem) {
482  qflux = k4 * phi_avg;
483  mlocal = (thermal.mlocal_ap) * thermal.massfraction + (thermal.mlocal_htpb) * (1.0 - thermal.massfraction);
484  mdota = fabs(mdot(i, j, k));
485  mbase = tanh(4.0 * mdota / (mlocal));
486  K = (thermal.modeling_ap * thermal.k_ap * thermal.massfraction + thermal.modeling_htpb * thermal.k_htpb * (1.0 - thermal.massfraction)) * phi_avg + thermal.disperssion1 * (1. - phi_avg);
487  heatflux(i, j, k) = (laser(i, j, k) * phi_avg + thermal.hc * mbase * qflux) / K;
488  }
489  else {
490  qflux = k1 * phi_avg +
491  k2 * (1.0 - phi_avg) +
492  (zeta_1 / zeta) * exp(k3 * phi_avg * (1.0 - phi_avg));
493  mlocal = (thermal.mlocal_ap) * phi_avg + (thermal.mlocal_htpb) * (1.0 - phi_avg) + thermal.mlocal_comb * phi_avg * (1.0 - phi_avg);
494  mdota = fabs(mdot(i, j, k));
495  mbase = tanh(4.0 * mdota / (mlocal));
496 
497  K = thermal.modeling_ap * thermal.k_ap * phi_avg + thermal.modeling_htpb * thermal.k_htpb * (1.0 - phi_avg);
498  heatflux(i, j, k) = (thermal.hc * mbase * qflux + laser(i, j, k)) / K;
499  }
500 
501  if (isnan(heatflux(i, j, k))) {
502  Util::Message(INFO, heatflux(i, j, k), "heat contains nan (i=", i, " j= ", j, ")");
503  Util::Abort(INFO);
504  }
505  });
506 
507  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
508  {
509  auto sten = Numeric::GetStencil(i, j, k, bx);
510  Set::Vector grad_eta = Numeric::Gradient(eta, i, j, k, 0, DX);
511  Set::Vector grad_temp = Numeric::Gradient(temp, i, j, k, 0, DX);
512  Set::Scalar lap_temp = Numeric::Laplacian(temp, i, j, k, 0, DX);
513  Set::Scalar grad_eta_mag = grad_eta.lpNorm<2>();
514  Set::Vector grad_alpha = Numeric::Gradient(alpha, i, j, k, 0, DX, sten);
515  Set::Scalar dTdt = 0.0;
516  dTdt += grad_eta.dot(grad_temp * alpha(i, j, k));
517  dTdt += grad_alpha.dot(eta(i, j, k) * grad_temp);
518  dTdt += eta(i, j, k) * alpha(i, j, k) * lap_temp;
519  dTdt += alpha(i, j, k) * heatflux(i, j, k) * grad_eta_mag;
520  Set::Scalar Tsolid;
521  Tsolid = dTdt + temps(i, j, k) * (etanew(i, j, k) - eta(i, j, k)) / dt;
522  tempsnew(i, j, k) = temps(i, j, k) + dt * Tsolid;
523  tempnew(i, j, k) = etanew(i, j, k) * tempsnew(i, j, k) + (1.0 - etanew(i, j, k)) * thermal.T_fluid;
524  if (isnan(tempsnew(i, j, k)) || isnan(temps(i, j, k))) {
525  Util::Message(INFO, tempsnew(i, j, k), "tempsnew contains nan (i=", i, " j= ", j, ")");
526  Util::Message(INFO, temps(i, j, k), "temps contains nan (i=", i, " j= ", j, ")");
527  Util::Abort(INFO);
528  }
529 
530  });
531 
532  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
533  {
534  Set::Scalar phi_avg = Numeric::Interpolate::NodeToCellAverage(phi, i, j, k, 0);
535  //phi_avg = phicell(i,j,k);
536 
537  Set::Scalar L;
538  if (pressure.arrhenius.mob_ap == 1) L = thermal.m_ap * pressure.P * exp(-thermal.E_ap / tempnew(i, j, k)) * phi_avg;
539  else L = thermal.m_ap * exp(-thermal.E_ap / tempnew(i, j, k)) * phi_avg;
540  L += thermal.m_htpb * exp(-thermal.E_htpb / tempnew(i, j, k)) * (1.0 - phi_avg);
541  //L += thermal.m_comb * (0.5 * tempnew(i, j, k) / thermal.bound) * phi(i, j, k) * (1.0 - phi(i, j, k));
542  if (tempnew(i, j, k) <= thermal.bound) mob(i, j, k) = 0;
543  else mob(i, j, k) = L;
544  //mob(i,j,k) = L;
545  if (isnan(mob(i, j, k))) {
546  Util::Message(INFO, mob(i, j, k), "mob contains nan (i=", i, " j= ", j, ")");
547  Util::Abort(INFO);
548  }
549  });
550  } // MFi For loop
551 
552  } // thermal IF
553  else
554  {
555  std::swap(eta_old_mf[lev], eta_mf[lev]);
557  0.0,
558  -5.0 * pf.w1 + 16.0 * pf.w12 - 11.0 * pf.w0,
559  14.0 * pf.w1 - 32.0 * pf.w12 + 18.0 * pf.w0,
560  -8.0 * pf.w1 + 16.0 * pf.w12 - 8.0 * pf.w0
561  );
563 
564  for (amrex::MFIter mfi(*eta_mf[lev], true); mfi.isValid(); ++mfi)
565  {
566  const amrex::Box& bx = mfi.tilebox();
567  //Phase Fields
568  amrex::Array4<Set::Scalar> const& etanew = (*eta_mf[lev]).array(mfi);
569  amrex::Array4<Set::Scalar> const& eta = (*eta_old_mf[lev]).array(mfi);
570  amrex::Array4<Set::Scalar> const& phi = (*phi_mf[lev]).array(mfi);
571 
572  Set::Scalar fmod_ap = pressure.power.r_ap * pow(pressure.P, pressure.power.n_ap);
573  Set::Scalar fmod_htpb = pressure.power.r_htpb * pow(pressure.P, pressure.power.n_htpb);
574  Set::Scalar fmod_comb = pressure.power.r_comb * pow(pressure.P, pressure.power.n_comb);
575 
576  pressure.power.a_fit = -1.16582 * sin(pressure.P) - 0.681788 * cos(pressure.P) + 3.3563;
577  pressure.power.b_fit = -0.708225 * sin(pressure.P) + 0.548067 * cos(pressure.P) + 1.55985;
578  pressure.power.c_fit = -0.0130849 * sin(pressure.P) - 0.03597 * cos(pressure.P) + 0.00725694;
579 
580 
581  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
582  {
583  //Set::Scalar phi_avg = Numeric::Interpolate::NodeToCellAverage(phi, i, j, k, 0);
584  Set::Scalar L;
585  if (homogeneousSystem == 1) {
586  Set::Vector grad_eta = Numeric::Gradient(eta, i, j, k, 0, DX);
587  Set::Scalar angle = acos(grad_eta[0] / grad_eta.lpNorm<2>()) * 180 / 3.1415;
588  if (angle > 90) angle = angle - 90.0;
589  if (angle > 180) angle = angle - 180.0;
590  if (angle > 270) angle = angle - 270.0;
591  L = pressure.power.a_fit + pressure.power.b_fit * exp(-pressure.power.c_fit * angle);
592  }
593  else {
594  Set::Scalar fs_actual;
595  fs_actual = fmod_ap * phi(i, j, k)
596  + fmod_htpb * (1.0 - phi(i, j, k))
597  + 4.0 * fmod_comb * phi(i, j, k) * (1.0 - phi(i, j, k));
598  L = fs_actual / pf.gamma / (pf.w1 - pf.w0);
599  }
600  Set::Scalar eta_lap = Numeric::Laplacian(eta, i, j, k, 0, DX);
601  Set::Scalar df_deta = ((pf.lambda / pf.eps) * dw(eta(i, j, k)) - pf.eps * pf.kappa * eta_lap);
602  etanew(i, j, k) = eta(i, j, k) - L * dt * df_deta;
603 
604  if (etanew(i, j, k) > eta(i, j, k)) etanew(i, j, k) = eta(i, j, k);
605  });
606  } // MFI For Loop
607  } // thermal ELSE
608  }// First IF
609 } //Function
610 
611 
612 void Flame::TagCellsForRefinement(int lev, amrex::TagBoxArray& a_tags, Set::Scalar time, int ngrow)
613 {
614  BL_PROFILE("Integrator::Flame::TagCellsForRefinement");
615  Base::Mechanics<model_type>::TagCellsForRefinement(lev, a_tags, time, ngrow);
616 
617  const Set::Scalar* DX = geom[lev].CellSize();
618  Set::Scalar dr = sqrt(AMREX_D_TERM(DX[0] * DX[0], +DX[1] * DX[1], +DX[2] * DX[2]));
619 
620  // Eta criterion for refinement
621  for (amrex::MFIter mfi(*eta_mf[lev], true); mfi.isValid(); ++mfi)
622  {
623  const amrex::Box& bx = mfi.tilebox();
624  amrex::Array4<char> const& tags = a_tags.array(mfi);
625  amrex::Array4<const Set::Scalar> const& eta = (*eta_mf[lev]).array(mfi);
626 
627  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
628  {
629  Set::Vector gradeta = Numeric::Gradient(eta, i, j, k, 0, DX);
630  if (gradeta.lpNorm<2>() * dr * 2 > m_refinement_criterion && eta(i, j, k) >= t_refinement_restriction)
631  tags(i, j, k) = amrex::TagBox::SET;
632  });
633  }
634 
635  // Phi criterion for refinement
636  if (elastic.phirefinement) {
637  for (amrex::MFIter mfi(*eta_mf[lev], true); mfi.isValid(); ++mfi)
638  {
639  const amrex::Box& bx = mfi.tilebox();
640  amrex::Array4<char> const& tags = a_tags.array(mfi);
641  amrex::Array4<const Set::Scalar> const& phi = (*phi_mf[lev]).array(mfi);
642 
643  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
644  {
645  Set::Vector gradphi = Numeric::Gradient(phi, i, j, k, 0, DX);
646  if (gradphi.lpNorm<2>() * dr >= phi_refinement_criterion)
647  tags(i, j, k) = amrex::TagBox::SET;
648  });
649  }
650  }
651 
652 
653  // Thermal criterion for refinement
654  if (thermal.on) {
655  for (amrex::MFIter mfi(*temp_mf[lev], true); mfi.isValid(); ++mfi)
656  {
657  const amrex::Box& bx = mfi.tilebox();
658  amrex::Array4<char> const& tags = a_tags.array(mfi);
659  amrex::Array4<const Set::Scalar> const& temp = (*temp_mf[lev]).array(mfi);
660  amrex::Array4<const Set::Scalar> const& eta = (*eta_mf[lev]).array(mfi);
661  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
662  {
663  Set::Vector tempgrad = Numeric::Gradient(temp, i, j, k, 0, DX);
664  if (tempgrad.lpNorm<2>() * dr > t_refinement_criterion && eta(i, j, k) >= t_refinement_restriction)
665  tags(i, j, k) = amrex::TagBox::SET;
666  });
667  }
668  }
669 
670 
671 }
672 
673 void Flame::Regrid(int lev, Set::Scalar time)
674 {
675  BL_PROFILE("Integrator::Flame::Regrid");
676  //if (lev < finest_level) return;
677  //phi_mf[lev]->setVal(0.0);
678  ic_phi->Initialize(lev, phi_mf, time);
679  //ic_phicell->Initialize(lev, phi_mf, time);
680 }
681 
682 void Flame::Integrate(int amrlev, Set::Scalar /*time*/, int /*step*/,
683  const amrex::MFIter& mfi, const amrex::Box& box)
684 {
685  BL_PROFILE("Flame::Integrate");
686  const Set::Scalar* DX = geom[amrlev].CellSize();
687  Set::Scalar dv = AMREX_D_TERM(DX[0], *DX[1], *DX[2]);
688  amrex::Array4<amrex::Real> const& eta = (*eta_mf[amrlev]).array(mfi);
689  amrex::Array4<amrex::Real> const& mdot = (*mdot_mf[amrlev]).array(mfi);
690  if (variable_pressure) {
691  amrex::ParallelFor(box, [=] AMREX_GPU_DEVICE(int i, int j, int k)
692  {
693  volume += eta(i, j, k, 0) * dv;
694  Set::Vector grad = Numeric::Gradient(eta, i, j, k, 0, DX);
695  Set::Scalar normgrad = grad.lpNorm<2>();
696  Set::Scalar da = normgrad * dv;
697  area += da;
698 
699  Set::Vector mgrad = Numeric::Gradient(mdot, i, j, k, 0, DX);
700  Set::Scalar mnormgrad = mgrad.lpNorm<2>();
701  Set::Scalar dm = mnormgrad * dv;
702  massflux += dm;
703 
704  });
705  }
706  else {
707  amrex::ParallelFor(box, [=] AMREX_GPU_DEVICE(int i, int j, int k)
708  {
709  volume += eta(i, j, k, 0) * dv;
710  Set::Vector grad = Numeric::Gradient(eta, i, j, k, 0, DX);
711  Set::Scalar normgrad = grad.lpNorm<2>();
712  Set::Scalar da = normgrad * dv;
713  area += da;
714  });
715  }
716  // time dependent pressure data from experimenta -> p = 0.0954521220950523 * exp(15.289993148880678 * t)
717 }
718 } // namespace Integrator
719 
720 
Integrator::Flame::hc
Set::Scalar hc
Definition: Flame.H:137
Numeric::CellGradientOnNode
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:673
IC::BMP
Definition: BMP.H:20
Integrator::Flame::power
struct Integrator::Flame::@2::@6 power
Integrator::Flame::ic_phi
IC::IC * ic_phi
Definition: Flame.H:80
Numeric::Interpolate::CellToNodeAverage
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:1256
Integrator::Flame
Definition: Flame.H:30
Integrator::Flame::modeling_htpb
Set::Scalar modeling_htpb
Definition: Flame.H:138
IC::IC::Initialize
void Initialize(const int &a_lev, Set::Field< Set::Scalar > &a_field, Set::Scalar a_time=0.0)
Definition: IC.H:26
Util::RealFillBoundary
AMREX_FORCE_INLINE void RealFillBoundary(amrex::FabArray< amrex::BaseFab< T >> &a_mf, const amrex::Geometry &)
Definition: Util.H:307
Integrator::Flame::Flame
Flame()
Definition: Flame.cpp:18
Integrator::Base::Mechanics::TimeStepBegin
virtual void TimeStepBegin(Set::Scalar a_time, int a_step) override
Definition: Mechanics.H:160
Integrator::Mechanics
Definition: Mechanics.H:101
Integrator::Flame::q0
Set::Scalar q0
Definition: Flame.H:133
Integrator::Flame::massflux
Set::Scalar massflux
Definition: Flame.H:161
Integrator::Base::Mechanics< model_type >
Integrator::Integrator::box
std::vector< amrex::Box > box
Definition: Integrator.H:367
Integrator::Flame::ic_laser
IC::IC * ic_laser
Definition: Flame.H:81
Integrator::Flame::mlocal_htpb
Set::Scalar mlocal_htpb
Definition: Flame.H:139
Integrator::Flame::disperssion2
Set::Scalar disperssion2
Definition: Flame.H:143
Integrator::Flame::rho_ap
Set::Scalar rho_ap
Definition: Flame.H:130
Integrator::Flame::m_htpb
Set::Scalar m_htpb
Definition: Flame.H:135
Integrator::Flame::y_len
Set::Scalar y_len
Definition: Flame.H:164
Integrator::Flame::TagCellsForRefinement
void TagCellsForRefinement(int lev, amrex::TagBoxArray &tags, amrex::Real, int) override
Definition: Flame.cpp:612
Integrator::Flame::mlocal_ap
Set::Scalar mlocal_ap
Definition: Flame.H:139
Integrator::Flame::model_htpb
model_type model_htpb
Definition: Flame.H:151
PSRead.H
IC::PSRead
Definition: PSRead.H:21
Integrator::Flame::laser_mf
Set::Field< Set::Scalar > laser_mf
Definition: Flame.H:76
Integrator::Flame::Integrate
void Integrate(int amrlev, Set::Scalar time, int step, const amrex::MFIter &mfi, const amrex::Box &box) override
Definition: Flame.cpp:682
Integrator::Flame::rho_htpb
Set::Scalar rho_htpb
Definition: Flame.H:130
Integrator::Base::Mechanics::psi_mf
Set::Field< Set::Scalar > psi_mf
Definition: Mechanics.H:455
Integrator::Flame::zeta
Set::Scalar zeta
Definition: Flame.H:87
Integrator::Flame::disperssion3
Set::Scalar disperssion3
Definition: Flame.H:144
Integrator::Flame::bc_eta
BC::BC< Set::Scalar > * bc_eta
Definition: Flame.H:79
Numeric::Function::Polynomial< 3 >
Definition: Function.H:13
Integrator::Flame::model_ap
model_type model_ap
Definition: Flame.H:151
Integrator::Base::Mechanics< Model::Solid::Finite::NeoHookeanPredeformed >::model_mf
Set::Field< Model::Solid::Finite::NeoHookeanPredeformed > model_mf
Definition: Mechanics.H:453
Integrator::Flame::Advance
void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
Definition: Flame.cpp:376
Integrator::Flame::small
Set::Scalar small
Definition: Flame.H:89
Integrator::Flame::chamber_area
Set::Scalar chamber_area
Definition: Flame.H:159
Expression.H
Integrator::Flame::T_fluid
Set::Scalar T_fluid
Definition: Flame.H:140
Integrator::Flame::Tref
Set::Scalar Tref
Definition: Flame.H:150
Integrator::Flame::temps_old_mf
Set::Field< Set::Scalar > temps_old_mf
Definition: Flame.H:64
Laminate.H
Integrator::Flame::plot_field
bool plot_field
Definition: Flame.H:94
Integrator::Base::Mechanics::m_type
Type m_type
Definition: Mechanics.H:459
Model::Solid::Finite::NeoHookeanPredeformed
Definition: NeoHookeanPredeformed.H:13
Integrator::Flame::eps
Set::Scalar eps
Definition: Flame.H:99
Set::Vector
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition: Base.H:20
Model::Solid::Finite::NeoHookeanPredeformed::F0
Set::Matrix F0
Definition: NeoHookeanPredeformed.H:38
Integrator::Flame::P
Set::Scalar P
Definition: Flame.H:108
Constant.H
Integrator::Flame::zeta_0
Set::Scalar zeta_0
Definition: Flame.H:88
ParmParse.H
Numeric::Gradient
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:655
Constant.H
pp_query
#define pp_query(...)
Definition: ParmParse.H:104
Integrator::Flame::mob_mf
Set::Field< Set::Scalar > mob_mf
Definition: Flame.H:68
Integrator::Flame::ic_eta
IC::IC * ic_eta
Definition: Flame.H:91
Integrator::Flame::massfraction
Set::Scalar massfraction
Definition: Flame.H:141
Integrator::Base::Mechanics::Initialize
void Initialize(int lev) override
Use the #ic object to initialize::Temp.
Definition: Mechanics.H:140
Integrator::Flame::gamma
Set::Scalar gamma
Definition: Flame.H:98
Integrator::Flame::m_ap
Set::Scalar m_ap
Definition: Flame.H:135
Flame.H
Integrator::Flame::base_time
Set::Scalar base_time
Definition: Flame.H:90
Integrator::Flame::UpdateModel
void UpdateModel(int a_step, Set::Scalar a_time) override
Definition: Flame.cpp:283
Numeric::Interpolate::NodeToCellAverage
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:1294
Set::Scalar
amrex::Real Scalar
Definition: Base.H:19
Integrator::Integrator::RegisterNewFab
void RegisterNewFab(Set::Field< Set::Scalar > &new_fab, BC::BC< Set::Scalar > *new_bc, int ncomp, int nghost, std::string name, bool writeout)
Definition: Integrator.cpp:292
Integrator::Flame::E_ap
Set::Scalar E_ap
Definition: Flame.H:136
pp_queryclass
#define pp_queryclass(...)
Definition: ParmParse.H:105
Integrator::Flame::ghost_count
int ghost_count
Definition: Flame.H:92
Integrator::Flame::kappa
Set::Scalar kappa
Definition: Flame.H:101
Integrator::Flame::w0
Set::Scalar w0
Definition: Flame.H:102
Integrator::Flame::t_refinement_restriction
Set::Scalar t_refinement_restriction
Definition: Flame.H:86
IC::PNG
Definition: PNG.H:24
Mechanics.H
Integrator::Flame::area
Set::Scalar area
Definition: Flame.H:158
Integrator::Flame::m_refinement_criterion
Set::Scalar m_refinement_criterion
Definition: Flame.H:84
Integrator::Flame::phi_mf
Set::Field< Set::Scalar > phi_mf
Definition: Flame.H:71
Integrator::Flame::eta_old_mf
Set::Field< Set::Scalar > eta_old_mf
Definition: Flame.H:67
Integrator::Flame::lambda
Set::Scalar lambda
Definition: Flame.H:100
Integrator::Flame::elastic
struct Integrator::Flame::@4 elastic
Integrator::Flame::Regrid
void Regrid(int lev, Set::Scalar time) override
Definition: Flame.cpp:673
Numeric::Laplacian
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])
Definition: Stencil.H:530
Integrator::Flame::bound
Set::Scalar bound
Definition: Flame.H:134
Integrator::Flame::mlocal_comb
Set::Scalar mlocal_comb
Definition: Flame.H:139
BC::Nothing
Definition: Nothing.H:11
Integrator::Flame::ic_temp
IC::IC * ic_temp
Definition: Flame.H:145
Integrator::Flame::variable_pressure
int variable_pressure
Definition: Flame.H:95
Integrator::Flame::TimeStepComplete
void TimeStepComplete(Set::Scalar a_time, int a_iter) override
Definition: Flame.cpp:364
Integrator::Flame::w12
Set::Scalar w12
Definition: Flame.H:102
Integrator::Flame::w1
Set::Scalar w1
Definition: Flame.H:102
Integrator::Flame::chamber_pressure
Set::Scalar chamber_pressure
Definition: Flame.H:160
Integrator::Flame::pf
struct Integrator::Flame::@1 pf
pp_query_default
#define pp_query_default(...)
Definition: ParmParse.H:100
Integrator::Flame::cp_ap
Set::Scalar cp_ap
Definition: Flame.H:132
Integrator::Flame::eta_mf
Set::Field< Set::Scalar > eta_mf
Definition: Flame.H:66
Util::Abort
void Abort(const char *msg)
Definition: Util.cpp:165
Integrator
Collection of numerical integrator objects.
Definition: AllenCahn.H:39
IC::Constant
Definition: Constant.H:18
Integrator::Flame::bc_temp
BC::BC< Set::Scalar > * bc_temp
Definition: Flame.H:78
Integrator::Flame::temp_mf
Set::Field< Set::Scalar > temp_mf
Definition: Flame.H:61
BMP.H
BC::Constant
Definition: Constant.H:89
Integrator::Flame::TimeStepBegin
void TimeStepBegin(Set::Scalar a_time, int a_iter) override
Definition: Flame.cpp:354
IC::Expression
Definition: Expression.H:44
Integrator::Flame::traction
Set::Scalar traction
Definition: Flame.H:152
Integrator::Integrator::RegisterNodalFab
void RegisterNodalFab(Set::Field< Set::Scalar > &new_fab, int ncomp, int nghost, std::string name, bool writeout)
Definition: Integrator.cpp:310
Integrator::Flame::x_len
Set::Scalar x_len
Definition: Flame.H:163
Integrator::Flame::alpha_mf
Set::Field< Set::Scalar > alpha_mf
Definition: Flame.H:73
Integrator::Flame::t_refinement_criterion
Set::Scalar t_refinement_criterion
Definition: Flame.H:85
Integrator::Flame::mdot_mf
Set::Field< Set::Scalar > mdot_mf
Definition: Flame.H:69
Integrator::Base::Mechanics::TagCellsForRefinement
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar, int) override
Definition: Mechanics.H:429
pp_query_validate
#define pp_query_validate(...)
Definition: ParmParse.H:101
Integrator::Base::Mechanics< Model::Solid::Finite::NeoHookeanPredeformed >::rho
Set::Scalar rho
Definition: Mechanics.H:472
pp_query_required
#define pp_query_required(...)
Definition: ParmParse.H:99
IO::ParmParse
Definition: ParmParse.H:110
Integrator::Flame::heatflux_mf
Set::Field< Set::Scalar > heatflux_mf
Definition: Flame.H:74
Integrator::Base::Mechanics::Advance
void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
Definition: Mechanics.H:238
Integrator::Flame::cp_htpb
Set::Scalar cp_htpb
Definition: Flame.H:132
Integrator::Flame::Parse
static void Parse(Flame &value, IO::ParmParse &pp)
Definition: Flame.cpp:27
Integrator::Base::Mechanics::plot_psi
bool plot_psi
Definition: Mechanics.H:500
Integrator::Flame::on
bool on
Definition: Flame.H:129
Integrator::Flame::thermal
struct Integrator::Flame::@3 thermal
Integrator::Flame::k_htpb
Set::Scalar k_htpb
Definition: Flame.H:131
Integrator::Flame::temps_mf
Set::Field< Set::Scalar > temps_mf
Definition: Flame.H:63
Integrator::Flame::arrhenius
struct Integrator::Flame::@2::@5 arrhenius
Integrator::Flame::pressure
struct Integrator::Flame::@2 pressure
Integrator::Flame::temp_old_mf
Set::Field< Set::Scalar > temp_old_mf
Definition: Flame.H:62
Integrator::Flame::E_htpb
Set::Scalar E_htpb
Definition: Flame.H:136
IC::Laminate
Definition: Laminate.H:16
Function.H
INFO
#define INFO
Definition: Util.H:20
Integrator::Flame::modeling_ap
Set::Scalar modeling_ap
Definition: Flame.H:138
Integrator::Integrator::dt
amrex::Vector< amrex::Real > dt
Timesteps for each level of refinement.
Definition: Integrator.H:303
Integrator::Flame::phi_refinement_criterion
Set::Scalar phi_refinement_criterion
Definition: Flame.H:83
Numeric::Function::Polynomial< 4 >::D
Polynomial< 3 > D() const
Definition: Function.H:42
Integrator::Flame::Initialize
void Initialize(int lev) override
Definition: Flame.cpp:240
Integrator::Integrator::RegisterIntegratedVariable
void RegisterIntegratedVariable(Set::Scalar *integrated_variable, std::string name, bool extensive=true)
Definition: Integrator.cpp:320
Numeric::GetStencil
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:20
Integrator::Base::Mechanics< Model::Solid::Finite::NeoHookeanPredeformed >::rhs_mf
Set::Field< Set::Vector > rhs_mf
Definition: Mechanics.H:462
Integrator::Base::Mechanics::psi_on
bool psi_on
Definition: Mechanics.H:456
Integrator::Flame::homogeneousSystem
int homogeneousSystem
Definition: Flame.H:93
Integrator::Flame::bc_psi
BC::BC< Set::Scalar > * bc_psi
Definition: Flame.H:166
Util::Message
void Message(std::string file, std::string func, int line, Args const &... args)
Definition: Util.H:133
Integrator::Flame::disperssion1
Set::Scalar disperssion1
Definition: Flame.H:142
Integrator::Flame::phirefinement
int phirefinement
Definition: Flame.H:153
Numeric::Function::Polynomial< 4 >
Definition: Function.H:30
Integrator::Flame::volume
Set::Scalar volume
Definition: Flame.H:157
Integrator::Flame::k_ap
Set::Scalar k_ap
Definition: Flame.H:131