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