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 :
18 : #include <cmath>
19 :
20 : namespace Integrator
21 : {
22 :
23 3 : Flame::Flame() : 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", false);
118 :
119 : // phase field initial condition
120 6 : pp.select<IC::Laminate,IC::Constant,IC::Expression,IC::BMP,IC::PNG>("pf.eta.ic",value.ic_eta,value.geom);
121 :
122 :
123 :
124 : // Select reduced order model to capture heat feedback
125 : pp.select< Model::Propellant::PowerLaw,
126 : Model::Propellant::FullFeedback,
127 : Model::Propellant::Homogenize>
128 6 : ("propellant",value.propellant);
129 :
130 :
131 : // Whether to use the Thermal Transport Model
132 6 : pp_query_default("thermal.on", value.thermal.on, false);
133 :
134 : // Reference temperature
135 : // Used to set all other reference temperatures by default.
136 12 : pp_query_default("thermal.Tref", value.thermal.Tref, "300.0_K",Unit::Temperature());
137 :
138 3 : if (value.thermal.on) {
139 :
140 : // Used to change heat flux units
141 4 : pp_query_default("thermal.hc", value.thermal.hc, "1.0", Unit::Power()/Unit::Area());
142 :
143 : // Effective fluid temperature
144 2 : pp_query_default("thermal.Tfluid", value.thermal.Tfluid, value.thermal.Tref);
145 :
146 : //Temperature boundary condition
147 2 : pp.select_default<BC::Constant>("thermal.temp.bc", value.bc_temp, 1, Unit::Temperature());
148 :
149 3 : value.RegisterNewFab(value.temp_mf, value.bc_temp, 1, 3, "temp", true);
150 3 : value.RegisterNewFab(value.temp_old_mf, value.bc_temp, 1, 3, "temp_old", false);
151 3 : value.RegisterNewFab(value.temps_mf, value.bc_temp, 1, 0, "temps", false);
152 :
153 3 : value.RegisterNewFab(value.mdot_mf, value.bc_temp, 1, 0, "mdot", value.plot_field);
154 3 : value.RegisterNewFab(value.alpha_mf, value.bc_temp, 1, 0, "alpha", value.plot_field);
155 3 : value.RegisterNewFab(value.heatflux_mf, value.bc_temp, 1, 0, "heatflux", value.plot_field);
156 3 : value.RegisterNewFab(value.laser_mf, value.bc_temp, 1, 0, "laser", value.plot_field);
157 :
158 2 : value.RegisterIntegratedVariable(&value.chamber.volume, "volume");
159 2 : value.RegisterIntegratedVariable(&value.chamber.area, "area");
160 2 : value.RegisterIntegratedVariable(&value.chamber.massflux, "mass_flux");
161 :
162 : // laser initial condition
163 : pp.select_default< IC::Constant,
164 : IC::Expression >
165 2 : ("laser.ic",value.ic_laser, value.geom, Unit::Power()/Unit::Area());
166 :
167 : // thermal initial condition
168 : pp.select_default< IC::Constant,
169 : IC::Expression,
170 : IC::BMP,
171 : IC::PNG >
172 3 : ("temp.ic",value.thermal.ic_temp,value.geom, Unit::Temperature());
173 : }
174 :
175 :
176 : // Constant pressure value
177 12 : pp_query_default("chamber.pressure", value.chamber.pressure, "1.0_MPa", Unit::Pressure());
178 :
179 : // Whether to compute the pressure evolution
180 6 : pp_query_default("variable_pressure", value.variable_pressure, false);
181 :
182 : // Refinement criterion for eta field
183 12 : pp_query_default( "amr.refinement_criterion", value.m_refinement_criterion, "0.001",
184 : Unit::Less());
185 :
186 : // Refinement criterion for temperature field
187 12 : pp.query_default( "amr.refinement_criterion_temp", value.t_refinement_criterion, "0.001_K",
188 : Unit::Temperature());
189 :
190 : // Eta value to restrict the refinament for the temperature field
191 12 : pp.query_default( "amr.refinament_restriction", value.t_refinement_restriction, "0.1",
192 : Unit::Less());
193 :
194 : // Refinement criterion for phi field [infinity]
195 6 : pp_query_default("amr.phi_refinement_criterion", value.phi_refinement_criterion, 1.0e100);
196 :
197 : // Minimum allowable threshold for $\eta$
198 6 : pp_query_default("small", value.small, 1.0e-8);
199 :
200 : // Initial condition for $\phi$ field.
201 : pp.select_default<IC::PSRead,IC::Laminate,IC::Expression,IC::Constant,IC::BMP,IC::PNG>
202 6 : ("phi.ic",value.ic_phi,value.geom);
203 :
204 9 : value.RegisterNodalFab(value.phi_mf, 1, 2, "phi", true);
205 :
206 : // Whether to use Neo-hookean Elastic model
207 6 : pp_query_default("elastic.on", value.elastic.on, 0);
208 :
209 : // Body force
210 6 : pp_query_default("elastic.traction", value.elastic.traction, 0.0);
211 :
212 : // Phi refinement criteria
213 6 : pp_query_default("elastic.phirefinement", value.elastic.phirefinement, 1);
214 :
215 6 : pp.queryclass<Base::Mechanics<model_type>>("elastic",value);
216 :
217 3 : if (value.m_type != Type::Disable)
218 : {
219 : // Reference temperature for thermal expansion
220 : // (temperature at which the material is strain-free)
221 0 : pp_query_default("Telastic", value.elastic.Telastic, value.thermal.Tref);
222 : // elastic model of AP
223 0 : pp.queryclass<Model::Solid::Finite::NeoHookeanPredeformed>("model_ap", value.elastic.model_ap);
224 : // elastic model of HTPB
225 0 : pp.queryclass<Model::Solid::Finite::NeoHookeanPredeformed>("model_htpb", value.elastic.model_htpb);
226 :
227 : // Use our current eta field as the psi field for the solver
228 0 : value.psi_on = false;
229 0 : value.solver.setPsi(value.eta_mf);
230 : }
231 :
232 : bool allow_unused;
233 : // Set this to true to allow unused inputs without error.
234 : // (Not recommended.)
235 3 : pp.query_default("allow_unused",allow_unused,false);
236 3 : if (!allow_unused && pp.AnyUnusedInputs())
237 : {
238 0 : Util::Warning(INFO,"The following inputs were specified but not used:");
239 0 : pp.AllUnusedInputs();
240 0 : Util::Exception(INFO,"Aborting. Specify 'allow_unused=True` to ignore this error.");
241 : }
242 3 : }
243 :
244 14 : void Flame::Initialize(int lev)
245 : {
246 : BL_PROFILE("Integrator::Flame::Initialize");
247 14 : Base::Mechanics<model_type>::Initialize(lev);
248 :
249 14 : ic_eta->Initialize(lev, eta_mf);
250 14 : ic_eta->Initialize(lev, eta_old_mf);
251 14 : ic_phi->Initialize(lev, phi_mf);
252 : //ic_phicell->Initialize(lev, phicell_mf);
253 :
254 14 : if (elastic.on) {
255 0 : rhs_mf[lev]->setVal(Set::Vector::Zero());
256 : }
257 14 : if (thermal.on) {
258 8 : if (thermal.ic_temp)
259 : {
260 8 : thermal.ic_temp->Initialize(lev,temp_mf);
261 8 : thermal.ic_temp->Initialize(lev,temp_old_mf);
262 8 : thermal.ic_temp->Initialize(lev,temps_mf);
263 : }
264 : else
265 : {
266 0 : temp_mf[lev]->setVal(thermal.Tref);
267 0 : temp_old_mf[lev]->setVal(thermal.Tref);
268 0 : temps_mf[lev]->setVal(thermal.Tref);
269 : }
270 8 : alpha_mf[lev]->setVal(0.0);
271 8 : mdot_mf[lev]->setVal(0.0);
272 8 : heatflux_mf[lev]->setVal(0.0);
273 8 : ic_laser->Initialize(lev, laser_mf);
274 : }
275 14 : if (variable_pressure) chamber.pressure = 1.0;
276 14 : }
277 :
278 0 : void Flame::UpdateModel(int /*a_step*/, Set::Scalar /*a_time*/)
279 : {
280 0 : if (m_type == Base::Mechanics<model_type>::Type::Disable) return;
281 :
282 0 : for (int lev = 0; lev <= finest_level; ++lev)
283 : {
284 0 : amrex::Box domain = this->geom[lev].Domain();
285 0 : domain.convert(amrex::IntVect::TheNodeVector());
286 0 : const Set::Scalar* DX = geom[lev].CellSize();
287 :
288 0 : phi_mf[lev]->FillBoundary();
289 0 : eta_mf[lev]->FillBoundary();
290 0 : temp_mf[lev]->FillBoundary();
291 :
292 0 : for (MFIter mfi(*model_mf[lev], false); mfi.isValid(); ++mfi)
293 : {
294 0 : amrex::Box smallbox = mfi.nodaltilebox();
295 0 : amrex::Box bx = mfi.grownnodaltilebox() & domain;
296 0 : Set::Patch<model_type> model = model_mf.Patch(lev,mfi);
297 0 : Set::Patch<const Set::Scalar> phi = phi_mf.Patch(lev,mfi);
298 0 : Set::Patch<const Set::Scalar> eta = eta_mf.Patch(lev,mfi);
299 0 : Set::Patch<Set::Vector> rhs = rhs_mf.Patch(lev,mfi);
300 :
301 0 : if (elastic.on)
302 : {
303 0 : Set::Patch <const Set::Scalar> temp = temp_mf.Patch(lev,mfi);
304 0 : amrex::ParallelFor(smallbox, [=] AMREX_GPU_DEVICE(int i, int j, int k)
305 : {
306 0 : Set::Vector grad_eta = Numeric::CellGradientOnNode(eta, i, j, k, 0, DX);
307 0 : rhs(i, j, k) = elastic.traction * grad_eta;
308 0 : });
309 0 : amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
310 : {
311 0 : Set::Scalar phi_avg = phi(i, j, k, 0);
312 0 : Set::Scalar temp_avg = Numeric::Interpolate::CellToNodeAverage(temp, i, j, k, 0);
313 0 : model_type model_ap = elastic.model_ap;
314 0 : model_ap.F0 -= Set::Matrix::Identity();
315 0 : model_ap.F0 *= (temp_avg - elastic.Telastic);
316 0 : model_ap.F0 += Set::Matrix::Identity();
317 0 : model_type model_htpb = elastic.model_htpb;
318 0 : model_htpb.F0 -= Set::Matrix::Identity();
319 0 : model_htpb.F0 *= (temp_avg - elastic.Telastic);
320 0 : model_htpb.F0 += Set::Matrix::Identity();
321 :
322 0 : model(i, j, k) = model_ap * phi_avg + model_htpb * (1. - phi_avg);
323 0 : });
324 : }
325 : else
326 : {
327 0 : amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
328 : {
329 0 : Set::Scalar phi_avg = Numeric::Interpolate::CellToNodeAverage(phi, i, j, k, 0);
330 : //phi_avg = phi(i,j,k,0);
331 0 : model_type model_ap = elastic.model_ap;
332 0 : model_ap.F0 *= Set::Matrix::Zero();
333 0 : model_type model_htpb = elastic.model_htpb;
334 0 : model_htpb.F0 *= Set::Matrix::Zero();
335 0 : model(i, j, k) = model_ap * phi_avg + model_htpb * (1. - phi_avg);
336 0 : });
337 : }
338 0 : }
339 0 : Util::RealFillBoundary(*model_mf[lev], geom[lev]);
340 :
341 : }
342 : }
343 :
344 818 : void Flame::TimeStepBegin(Set::Scalar a_time, int a_iter)
345 : {
346 : BL_PROFILE("Integrator::Flame::TimeStepBegin");
347 818 : Base::Mechanics<model_type>::TimeStepBegin(a_time, a_iter);
348 818 : if (thermal.on) {
349 90 : for (int lev = 0; lev <= finest_level; ++lev)
350 80 : ic_laser->Initialize(lev, laser_mf, a_time);
351 : }
352 818 : }
353 :
354 818 : void Flame::TimeStepComplete(Set::Scalar /*a_time*/, int /*a_iter*/)
355 : {
356 : BL_PROFILE("Integrator::Flame::TimeStepComplete");
357 818 : if (variable_pressure) {
358 : //Set::Scalar x_len = geom[0].ProbDomain().length(0);
359 : //Set::Scalar y_len = geom[0].ProbDomain().length(1);
360 : // Set::Scalar domain_area = x_len * y_len;
361 0 : Util::Message(INFO, "Mass = ", chamber.massflux);
362 0 : Util::Message(INFO, "Pressure = ", chamber.pressure);
363 : }
364 818 : }
365 :
366 8206 : void Flame::Advance(int lev, Set::Scalar time, Set::Scalar dt)
367 : {
368 : BL_PROFILE("Integrador::Flame::Advance");
369 8206 : Base::Mechanics<model_type>::Advance(lev, time, dt);
370 8206 : const Set::Scalar* DX = geom[lev].CellSize();
371 :
372 8206 : std::swap(eta_old_mf[lev], eta_mf[lev]);
373 :
374 : //
375 : // Chamber pressure update
376 : //
377 8206 : if (variable_pressure) {
378 0 : chamber.pressure = exp(0.00075 * chamber.massflux);
379 0 : if (chamber.pressure > 10.0) {
380 0 : chamber.pressure = 10.0;
381 : }
382 0 : else if (chamber.pressure <= 0.99) {
383 0 : chamber.pressure = 0.99;
384 : }
385 0 : elastic.traction = chamber.pressure;
386 : }
387 :
388 :
389 : //
390 : // Multi-well chemical potential
391 : //
392 : Numeric::Function::Polynomial<4> w( pf.w0,
393 : 0.0,
394 8206 : -5.0 * pf.w1 + 16.0 * pf.w12 - 11.0 * pf.w0,
395 8206 : 14.0 * pf.w1 - 32.0 * pf.w12 + 18.0 * pf.w0,
396 8206 : -8.0 * pf.w1 + 16.0 * pf.w12 - 8.0 * pf.w0);
397 8206 : Numeric::Function::Polynomial<3> dw = w.D();
398 :
399 8206 : propellant.set_pressure(chamber.pressure);
400 :
401 30068 : for (amrex::MFIter mfi(*eta_mf[lev], true); mfi.isValid(); ++mfi)
402 : {
403 21862 : const amrex::Box& bx = mfi.tilebox();
404 : // Phase fields
405 21862 : Set::Patch<Set::Scalar> etanew = eta_mf.Patch(lev,mfi);
406 21862 : Set::Patch<const Set::Scalar> eta = eta_old_mf.Patch(lev,mfi);
407 21862 : Set::Patch<const Set::Scalar> phi = phi_mf.Patch(lev,mfi);
408 : // Heat transfer fields
409 21862 : Set::Patch<const Set::Scalar> temp = temp_mf.Patch(lev,mfi);
410 21862 : Set::Patch<Set::Scalar> alpha = alpha_mf.Patch(lev,mfi);
411 21862 : Set::Patch<Set::Scalar> laser = laser_mf.Patch(lev,mfi);
412 : // Diagnostic fields
413 21862 : Set::Patch<Set::Scalar> mdot = mdot_mf.Patch(lev,mfi);
414 21862 : Set::Patch<Set::Scalar> heatflux = heatflux_mf.Patch(lev,mfi);
415 :
416 :
417 21862 : amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
418 : {
419 : //
420 : // CALCULATE PHI-AVERAGED QUANTITIES
421 : //
422 59108224 : Set::Scalar phi_avg = Numeric::Interpolate::NodeToCellAverage(phi, i, j, k, 0);
423 60856064 : Set::Scalar T = thermal.on ? temp(i,j,k) : NAN;
424 :
425 59108224 : Set::Scalar K = propellant.get_K(phi_avg);
426 :
427 59108224 : Set::Scalar rho = propellant.get_rho(phi_avg);
428 :
429 59108224 : Set::Scalar cp = propellant.get_cp(phi_avg);
430 :
431 : //
432 : // CALCULATE MOBILITY
433 : //
434 59108224 : Set::Scalar L = propellant.get_L( phi_avg, T);
435 :
436 : //
437 : // EVOLVE PHASE FIELD (ETA)
438 : //
439 :
440 59108224 : Set::Scalar eta_lap = Numeric::Laplacian(eta, i, j, k, 0, DX);
441 118216448 : Set::Scalar df_deta = ((pf.lambda / pf.eps) * dw(eta(i, j, k)) - pf.eps * pf.kappa * eta_lap);
442 :
443 59108224 : if (df_deta < 0) {
444 : // Prevent eta from increasing/healing. A bug was found where if the diffuse thickness was too large compared to a void
445 : // (region of eta = 0), eta would heal/increase in a non-physcial way, this if statement stops that behavior
446 196236 : df_deta = 0.0;
447 : }
448 :
449 118216448 : etanew(i, j, k) = eta(i, j, k) - L * dt * df_deta;
450 118216448 : if (etanew(i, j, k) <= small) etanew(i, j, k) = small;
451 :
452 59108224 : if (thermal.on)
453 : {
454 : //
455 : // Calculate thermal diffisivity and store for later gradient
456 : //
457 :
458 1747840 : alpha(i, j, k) = K / rho / cp;
459 :
460 : //
461 : // CALCULATE MASS FLUX BASED ON EVOLVING ETA
462 : //
463 :
464 5243520 : mdot(i, j, k) = rho * fabs(eta(i, j, k) - etanew(i, j, k)) / dt;
465 :
466 : //
467 : // CALCULATE HEAT FLUX BASED ON THE CALCULATED MASS FLUX
468 : //
469 :
470 3495680 : Set::Scalar q0 = propellant.get_qdot(mdot(i,j,k), phi_avg);
471 5243520 : heatflux(i,j,k) = ( thermal.hc*q0 + laser(i,j,k) ) / K;
472 :
473 : }
474 :
475 59108224 : });
476 :
477 8206 : } // MFi For loop
478 :
479 :
480 : //
481 : // THERMAL TRANSPORT
482 : //
483 8206 : if (thermal.on)
484 : {
485 2550 : std::swap(temp_old_mf[lev], temp_mf[lev]);
486 :
487 13100 : for (amrex::MFIter mfi(*eta_mf[lev], true); mfi.isValid(); ++mfi)
488 : {
489 10550 : const amrex::Box& bx = mfi.tilebox();
490 :
491 10550 : Set::Patch<Set::Scalar> tempnew = temp_mf.Patch(lev,mfi);
492 10550 : Set::Patch<const Set::Scalar> temp = temp_old_mf.Patch(lev,mfi);
493 10550 : Set::Patch<const Set::Scalar> alpha = alpha_mf.Patch(lev,mfi);
494 :
495 10550 : Set::Patch<Set::Scalar> temps = temps_mf.Patch(lev,mfi);
496 :
497 :
498 : // Phase field
499 10550 : Set::Patch<Set::Scalar> etanew = (*eta_mf[lev]).array(mfi);
500 10550 : Set::Patch<const Set::Scalar> eta = (*eta_old_mf[lev]).array(mfi);
501 : // Diagnostic fields
502 10550 : Set::Patch<const Set::Scalar> heatflux = heatflux_mf.Patch(lev,mfi);
503 :
504 10550 : amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
505 : {
506 1747840 : auto sten = Numeric::GetStencil(i, j, k, bx);
507 1747840 : Set::Vector grad_eta = Numeric::Gradient(eta, i, j, k, 0, DX);
508 1747840 : Set::Vector grad_temp = Numeric::Gradient(temp, i, j, k, 0, DX);
509 1747840 : Set::Scalar lap_temp = Numeric::Laplacian(temp, i, j, k, 0, DX);
510 1747840 : Set::Scalar grad_eta_mag = grad_eta.lpNorm<2>();
511 1747840 : Set::Vector grad_alpha = Numeric::Gradient(alpha, i, j, k, 0, DX, sten);
512 1747840 : Set::Scalar dTdt = 0.0;
513 3495680 : dTdt += grad_eta.dot(grad_temp * alpha(i, j, k));
514 3495680 : dTdt += grad_alpha.dot(eta(i, j, k) * grad_temp);
515 3495680 : dTdt += eta(i, j, k) * alpha(i, j, k) * lap_temp;
516 3495680 : dTdt += alpha(i, j, k) * heatflux(i, j, k) * grad_eta_mag;
517 :
518 5243520 : Set::Scalar Tsolid = dTdt + temps(i, j, k) * (etanew(i, j, k) - eta(i, j, k)) / dt;
519 3495680 : temps(i, j, k) = temps(i, j, k) + dt * Tsolid;
520 6991360 : tempnew(i, j, k) = etanew(i, j, k) * temps(i, j, k) + (1.0 - etanew(i, j, k)) * thermal.Tfluid;
521 1747840 : });
522 2550 : }
523 : }
524 :
525 8206 : } //Function
526 :
527 :
528 190 : void Flame::TagCellsForRefinement(int lev, amrex::TagBoxArray& a_tags, Set::Scalar time, int ngrow)
529 : {
530 : BL_PROFILE("Integrator::Flame::TagCellsForRefinement");
531 190 : Base::Mechanics<model_type>::TagCellsForRefinement(lev, a_tags, time, ngrow);
532 :
533 190 : const Set::Scalar* DX = geom[lev].CellSize();
534 190 : Set::Scalar dr = sqrt(AMREX_D_TERM(DX[0] * DX[0], +DX[1] * DX[1], +DX[2] * DX[2]));
535 :
536 : // Eta criterion for refinement
537 602 : for (amrex::MFIter mfi(*eta_mf[lev], true); mfi.isValid(); ++mfi)
538 : {
539 412 : const amrex::Box& bx = mfi.tilebox();
540 412 : amrex::Array4<char> const& tags = a_tags.array(mfi);
541 412 : Set::Patch<const Set::Scalar> eta = eta_mf.Patch(lev,mfi);
542 :
543 412 : amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
544 : {
545 519920 : Set::Vector gradeta = Numeric::Gradient(eta, i, j, k, 0, DX);
546 541390 : if (gradeta.lpNorm<2>() * dr * 2 > m_refinement_criterion && eta(i, j, k) >= t_refinement_restriction)
547 41936 : tags(i, j, k) = amrex::TagBox::SET;
548 519920 : });
549 190 : }
550 :
551 : // Phi criterion for refinement
552 190 : if (elastic.phirefinement) {
553 602 : for (amrex::MFIter mfi(*eta_mf[lev], true); mfi.isValid(); ++mfi)
554 : {
555 412 : const amrex::Box& bx = mfi.tilebox();
556 412 : amrex::Array4<char> const& tags = a_tags.array(mfi);
557 412 : Set::Patch<const Set::Scalar> phi = phi_mf.Patch(lev,mfi);
558 :
559 412 : amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
560 : {
561 519920 : Set::Vector gradphi = Numeric::Gradient(phi, i, j, k, 0, DX);
562 519920 : if (gradphi.lpNorm<2>() * dr >= phi_refinement_criterion)
563 0 : tags(i, j, k) = amrex::TagBox::SET;
564 519920 : });
565 190 : }
566 : }
567 :
568 :
569 : // Thermal criterion for refinement
570 190 : if (thermal.on) {
571 507 : for (amrex::MFIter mfi(*temp_mf[lev], true); mfi.isValid(); ++mfi)
572 : {
573 345 : const amrex::Box& bx = mfi.tilebox();
574 345 : amrex::Array4<char> const& tags = a_tags.array(mfi);
575 345 : Set::Patch<const Set::Scalar> temp = temp_mf.Patch(lev,mfi);
576 345 : Set::Patch<const Set::Scalar> eta = eta_mf.Patch(lev,mfi);
577 345 : amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
578 : {
579 50928 : Set::Vector tempgrad = Numeric::Gradient(temp, i, j, k, 0, DX);
580 50928 : if (tempgrad.lpNorm<2>() * dr > t_refinement_criterion && eta(i, j, k) >= t_refinement_restriction)
581 0 : tags(i, j, k) = amrex::TagBox::SET;
582 50928 : });
583 162 : }
584 : }
585 190 : }
586 :
587 2 : void Flame::Regrid(int lev, Set::Scalar time)
588 : {
589 : BL_PROFILE("Integrator::Flame::Regrid");
590 : //if (lev < finest_level) return;
591 : //phi_mf[lev]->setVal(0.0);
592 2 : ic_phi->Initialize(lev, phi_mf, time);
593 : //ic_phicell->Initialize(lev, phi_mf, time);
594 2 : }
595 :
596 220 : void Flame::Integrate(int amrlev, Set::Scalar /*time*/, int /*step*/,
597 : const amrex::MFIter& mfi, const amrex::Box& box)
598 : {
599 : BL_PROFILE("Flame::Integrate");
600 220 : const Set::Scalar* DX = geom[amrlev].CellSize();
601 220 : Set::Scalar dv = AMREX_D_TERM(DX[0], *DX[1], *DX[2]);
602 220 : Set::Patch<const Set::Scalar> eta = eta_mf.Patch(amrlev,mfi);
603 220 : Set::Patch<const Set::Scalar> mdot = mdot_mf.Patch(amrlev,mfi);
604 220 : if (variable_pressure) {
605 0 : amrex::ParallelFor(box, [=] AMREX_GPU_DEVICE(int i, int j, int k)
606 : {
607 0 : chamber.volume += eta(i, j, k, 0) * dv;
608 0 : Set::Vector grad = Numeric::Gradient(eta, i, j, k, 0, DX);
609 0 : Set::Scalar normgrad = grad.lpNorm<2>();
610 0 : Set::Scalar da = normgrad * dv;
611 0 : chamber.area += da;
612 :
613 0 : Set::Vector mgrad = Numeric::Gradient(mdot, i, j, k, 0, DX);
614 0 : Set::Scalar mnormgrad = mgrad.lpNorm<2>();
615 0 : Set::Scalar dm = mnormgrad * dv;
616 0 : chamber.massflux += dm;
617 :
618 0 : });
619 : }
620 : else {
621 220 : amrex::ParallelFor(box, [=] AMREX_GPU_DEVICE(int i, int j, int k)
622 : {
623 15560 : chamber.volume += eta(i, j, k, 0) * dv;
624 15560 : Set::Vector grad = Numeric::Gradient(eta, i, j, k, 0, DX);
625 15560 : Set::Scalar normgrad = grad.lpNorm<2>();
626 15560 : Set::Scalar da = normgrad * dv;
627 15560 : chamber.area += da;
628 15560 : });
629 : }
630 : // time dependent pressure data from experimenta -> p = 0.0954521220950523 * exp(15.289993148880678 * t)
631 220 : }
632 : } // namespace Integrator
633 :
634 :
|