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