Alamo
Mechanics.H
Go to the documentation of this file.
1#ifndef INTEGRATOR_BASE_MECHANICS_H
2#define INTEGRATOR_BASE_MECHANICS_H
3
4#include "AMReX.H"
10#include "Numeric/Stencil.H"
11#include "Model/Solid/Solid.H"
14#include "Operator/Operator.H"
15#include "IC/Constant.H"
16#include "IC/Expression.H"
17#include "IC/Trig.H"
18
19namespace Integrator
20{
21namespace Base
22{
23template<class MODEL>
24class Mechanics: virtual public Integrator
25{
26public:
27
29
31 {}
32
34 {
35 delete ic_rhs;
36 delete velocity_ic;
37 delete bc;
38 }
39
40 // The mechanics integrator manages the solution of an elastic
41 // solve using the MLMG solver.
42 static void Parse(Mechanics& value, IO::ParmParse& pp)
43 {
44 BL_PROFILE("Integrator::Base::Mechanics::Parse");
45 if (pp.contains("type"))
46 {
47 std::string type_str;
48 // Type of mecahnics to use.
49 // Static: do full implicit solve.
50 // Dynamic: evolve dynamic equations with explicit dynamics
51 // Disable: do nothing.
52 pp_query_validate("type", type_str, {"disable","static","dynamic"});
53 if (type_str == "static") value.m_type = Type::Static;
54 else if (type_str == "dynamic") value.m_type = Type::Dynamic;
55 else if (type_str == "disable") value.m_type = Type::Disable;
56 else Util::Abort(INFO, "Invalid type ", type_str, " specified");
57 }
58 if (value.m_type == Type::Disable) return;
59
60 // Treat mechanics fields as changing in time. [false]
61 // You should use this if you care about other physics driven by
62 // the output of this integrator.
63 pp_query_default("time_evolving", value.m_time_evolving,false);
64
65 pp_query_default("plot_disp", value.plot_disp, true); // Include displacement field in output
66 pp_query_default("plot_rhs", value.plot_rhs, true); // Include right-hand side in output
67 pp_query_default("plot_psi", value.plot_psi, true); // Include :math:`\psi` field in output
68 pp_query_default("plot_stress", value.plot_stress, true); // Include stress in output
69 pp_query_default("plot_strain", value.plot_strain, true); // Include strain in output
70
71 value.RegisterGeneralFab(value.disp_mf, 1, 2, value.plot_disp, "disp", value.m_time_evolving);
72 value.RegisterGeneralFab(value.rhs_mf, 1, 2, value.plot_rhs, "rhs", value.m_time_evolving);
73 value.RegisterGeneralFab(value.stress_mf, 1, 2, value.plot_stress, "stress", value.m_time_evolving);
74 value.RegisterGeneralFab(value.strain_mf, 1, 2, value.plot_strain, "strain", value.m_time_evolving);
75
76 if (value.m_type == Type::Static)
77 {
78 // Read parameters for :ref:`Solver::Nonlocal::Newton` solver
79 pp_queryclass("solver", value.solver);
80 }
81 if (value.m_type == Type::Dynamic)
82 {
83 value.RegisterGeneralFab(value.vel_mf, 1, 2, "vel",true);
84 value.RegisterGeneralFab(value.disp_old_mf, 1, 2, "dispold");
85 value.RegisterGeneralFab(value.vel_old_mf, 1, 2, "velold");
86 value.RegisterGeneralFab(value.ddw_mf, 1, 2);
87 pp_forbid("viscous.mu","replaced with viscous.mu_dashpot");
88 pp_query_default("viscous.mu_dashpot", value.mu_dashpot,0.0); // Dashpot damping (damps velocity)
89 pp_forbid("viscous.mu2","replaced with viscous.mu_newton");
90 pp_query_default("viscous.mu_newton", value.mu_newton,0.0); // Newtonian viscous damping (damps velocity gradient)
91
92 std::string velocity_ic_str;
93 pp_query_validate("velocity.ic.type",velocity_ic_str,{"none","expression"}); // Initializer for RHS
94 if (velocity_ic_str == "expression") value.velocity_ic = new IC::Expression(value.geom, pp,"velocity.ic.expression");
95 }
96
97 // Select the mechanical boundary conditions
99
100 pp_query_default("print_model", value.m_print_model, false); // Print out model variables (if enabled by model)
101 if (value.m_print_model) value.RegisterGeneralFab(value.model_mf, 1, 2, "model", value.m_time_evolving);
102 else value.RegisterGeneralFab(value.model_mf, 1, 2, value.m_time_evolving);
103
104 // This if is here because the code currently allows for no
105 // body force to be defined, in which case value.ic_rhs is set to
106 // a nullptr. This should eventually be fixed - this is a stopgap.
107 if (pp.contains("rhs.type"))
108 // initial condition for right hand side (body force)
109 pp.select<IC::Constant,IC::Expression,IC::Trig>("rhs",value.ic_rhs,value.geom);
110
111 // Timestep interval for elastic solves (default - solve every time)
112 pp_query_default("interval", value.m_interval, 0);
113
114 value.RegisterIntegratedVariable(&(value.disp_hi[0].data()[0]), "disp_xhi_x");
115 value.RegisterIntegratedVariable(&(value.disp_hi[0].data()[1]), "disp_xhi_y");
116 value.RegisterIntegratedVariable(&(value.disp_hi[1].data()[0]), "disp_yhi_x");
117 value.RegisterIntegratedVariable(&(value.disp_hi[1].data()[1]), "disp_yhi_y");
118 value.RegisterIntegratedVariable(&(value.trac_hi[0].data()[0]), "trac_xhi_x");
119 value.RegisterIntegratedVariable(&(value.trac_hi[0].data()[1]), "trac_xhi_y");
120 value.RegisterIntegratedVariable(&(value.trac_hi[1].data()[0]), "trac_yhi_x");
121 value.RegisterIntegratedVariable(&(value.trac_hi[1].data()[1]), "trac_yhi_y");
122
123 // Maximum multigrid coarsening level (default - none, maximum coarsening)
124 pp_query_default("max_coarsening_level", value.m_max_coarsening_level,-1);
125
126 // Whether to include residual output field
127 pp_query_default("print_residual", value.m_print_residual,false);
128 if (value.m_print_residual) value.RegisterGeneralFab(value.res_mf, 1, 2, "res", false);
129
130 // Whether to refine based on elastic solution
131 pp_query_default("elastic_ref_threshold", value.m_elastic_ref_threshold,0.01);
132
133 // Set this to true to zero out the displacement before each solve.
134 // (This is a temporary fix - we need to figure out why this is needed.)
135 pp_query_default("zero_out_displacement", value.m_zero_out_displacement,false);
136
137 // Time to start doing the elastic solve (by default, start immediately)
138 pp_query_default("tstart", value.tstart,-1.0);
139
140
141 // Relative tolerance in mechanics solve
142 pp.query_default("tol_rel",value.tol_rel,1E-8);
143
144 // Absolute tolerance in mechanics solve
145 pp.query_default("tol_abs",value.tol_abs,1E-8);
146 }
147
148protected:
149 /// \brief Use the #ic object to initialize#Temp
150 void Initialize(int lev) override
151 {
152 BL_PROFILE("Integrator::Base::Mechanics::Initialize");
154
155 disp_mf[lev]->setVal(Set::Vector::Zero());
156 //disp_old_mf[lev]->setVal(Set::Vector::Zero());
157
158 if (m_type == Type::Dynamic)
159 if (velocity_ic)
160 {
162 }
163
164 if (ic_rhs) ic_rhs->Initialize(lev, rhs_mf);
165 else rhs_mf[lev]->setVal(Set::Vector::Zero());
166 }
167
168 virtual void UpdateModel(int a_step, Set::Scalar a_time) = 0;
169
170 virtual void TimeStepBegin(Set::Scalar a_time, int a_step) override
171 {
172 BL_PROFILE("Integrator::Base::Mechanics::TimeStepBegin");
174
175 for (int lev = 0; lev <= finest_level; ++lev)
176 {
177 rhs_mf[lev]->setVal(Set::Vector::Zero());
178 if (ic_rhs) ic_rhs->Initialize(lev, rhs_mf);
179 }
180
181 UpdateModel(a_step, a_time);
182
183 bc->SetTime(a_time);
184 bc->Init(rhs_mf, geom);
185
187 if (a_time < tstart) return;
188 if (m_interval && a_step % m_interval) return;
189
190 amrex::LPInfo info;
191 if (m_max_coarsening_level >= 0)
192 info.setMaxCoarseningLevel(m_max_coarsening_level);
193 Operator::Elastic<MODEL::sym> elastic_op(Geom(0, finest_level), grids, DistributionMap(0, finest_level), info);
194 elastic_op.SetUniform(false);
195 elastic_op.SetHomogeneous(false);
196 elastic_op.SetBC(bc);
197 IO::ParmParse pp("elasticop");
198 pp_queryclass(elastic_op);
199
200 solver.Define(elastic_op);
202
203 for (int lev = 0; lev <= finest_level; ++lev)
204 if (m_zero_out_displacement) disp_mf[lev]->setVal(Set::Vector::Zero());
205
208 solver.Clear();
209
210 for (int lev = 0; lev <= disp_mf.finest_level; lev++)
211 {
212 amrex::Box domain = geom[lev].Domain();
213 domain.convert(amrex::IntVect::TheNodeVector());
214
215 const amrex::Real* DX = geom[lev].CellSize();
216 for (MFIter mfi(*disp_mf[lev], false); mfi.isValid(); ++mfi)
217 {
218 amrex::Box bx = mfi.nodaltilebox();
219 bx.grow(2);
220 bx = bx & domain;
221 amrex::Array4<MODEL> const& model = model_mf[lev]->array(mfi);
222 amrex::Array4<Set::Matrix> const& stress = stress_mf[lev]->array(mfi);
223 amrex::Array4<Set::Matrix> const& strain = strain_mf[lev]->array(mfi);
224 amrex::Array4<const Set::Vector> const& disp = disp_mf[lev]->array(mfi);
225
226
227 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
228 {
229 auto sten = Numeric::GetStencil(i, j, k, bx);
230 if (model(i, j, k).kinvar == Model::Solid::KinematicVariable::F)
231 {
232 Set::Matrix F = Set::Matrix::Identity() + Numeric::Gradient(disp, i, j, k, DX, sten);
233 stress(i, j, k) = model(i, j, k).DW(F);
234 strain(i, j, k) = F;
235 }
236 else
237 {
238 Set::Matrix gradu = Numeric::Gradient(disp, i, j, k, DX, sten);
239 stress(i, j, k) = model(i, j, k).DW(gradu);
240 strain(i, j, k) = 0.5 * (gradu + gradu.transpose());
241 }
242 });
243 }
244 stress_mf[lev]->setMultiGhost(true);
245 stress_mf[lev]->FillBoundaryAndSync(geom[lev].periodicity());
246 strain_mf[lev]->setMultiGhost(true);
247 strain_mf[lev]->FillBoundaryAndSync(geom[lev].periodicity());
248 }
249 }
250
251 void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
252 {
253 BL_PROFILE("Integrator::Base::Mechanics::Advance");
255 const amrex::Real* DX = geom[lev].CellSize();
256
257 amrex::Box domain = geom[lev].Domain();
258 domain.convert(amrex::IntVect::TheNodeVector());
259 const amrex::Dim3 lo = amrex::lbound(domain), hi = amrex::ubound(domain);
260
261 if (m_type == Type::Dynamic)
262 {
263
264 std::swap(*disp_mf[lev], *disp_old_mf[lev]);
265 std::swap(*vel_mf[lev], *vel_old_mf[lev]);
266
267
268 for (amrex::MFIter mfi(*disp_mf[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
269 {
273
274 Set::Patch<Set::Vector> unew = disp_mf.Patch(lev,mfi);
275 Set::Patch<Set::Vector> vnew = vel_mf.Patch(lev,mfi);
278 //Set::Patch<MATRIX4> ddw = ddw_mf.Patch(lev,mfi);
279 Set::Patch<MODEL> model = model_mf.Patch(lev,mfi);
280
281 amrex::Box bx = mfi.grownnodaltilebox() & domain;
282 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
283 {
284 auto sten = Numeric::GetStencil(i,j,k,bx);
285 eps(i, j, k) = Numeric::Gradient(u, i, j, k, DX,sten);
286 sig(i, j, k) = model(i, j, k).DW(eps(i, j, k));
287 //ddw(i,j,k) = model(i,j,k).DDW(eps(i,j,k));
288 });
289
290 bx = mfi.nodaltilebox() & domain;
291 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
292 {
293
294 bool AMREX_D_DECL(xmin = (i == lo.x), ymin = (j == lo.y), zmin = (k == lo.z));
295 bool AMREX_D_DECL(xmax = (i == hi.x), ymax = (j == hi.y), zmax = (k == hi.z));
296
297 if (AMREX_D_TERM(xmax || xmin, || ymax || ymin, || zmax || zmin))
298 {
299 auto sten = Numeric::GetStencil(i,j,k,domain);
300
301 auto bctype = bc->getType(i,j,k,domain);
302
303 for (int d = 0; d < AMREX_SPACEDIM; d++)
304 {
306 {
307 unew(i,j,k)(d) = b(i,j,k)(d);
308 }
310 {
311
312 Set::Vector N = Set::Normal(AMREX_D_DECL(xmin,ymin,zmin),
313 AMREX_D_DECL(xmax,ymax,zmax));
314
315 auto [phi, offdiag] = Numeric::GradientSplit(u,i,j,k,DX,sten);
316
317 Set::Matrix A = Set::Matrix::Zero();
318 Set::Vector rhs = b(i,j,k);
319 Set::Matrix DW_F0 = model(i,j,k).DW(Set::Matrix::Zero());
320 MATRIX4 ddw = model(i,j,k).DDW(eps(i,j,k));
321 //Util::Message(INFO,b(i,j,k).transpose());
322
323 for (int p = 0; p < AMREX_SPACEDIM; p++)
324 for (int q = 0; q < AMREX_SPACEDIM; q++)
325 {
326 for (int r = 0; r < AMREX_SPACEDIM; r++)
327 for (int s = 0; s < AMREX_SPACEDIM; s++)
328 {
329 A(p,r) += ddw(p,q,r,s) * phi(s) * N(q);
330
331 rhs(p) -= ddw(p,q,r,s) * eps(i,j,k)(r,s) * N(q);
332 }
333
334 rhs(p) -= DW_F0(p,q)*N(q);
335 }
336 Set::Vector delta_u = (A.inverse() * rhs);
337
338 unew(i,j,k)(d) = u(i,j,k)(d) + delta_u(d);
339 vnew(i,j,k)(d) = (unew(i,j,k)(d) - u(i,j,k)(d))/dt;
340 }
341 else
342 {
343 Util::Abort(INFO,"Elastic dynamics not supported for other BCs yet");
344 }
345 }
346 }
347 else
348 {
349 //Set::Matrix gradu = Numeric::Gradient(u,i,j,k,DX);
350 //Set::Matrix3 gradgradu = Numeric::Hessian(u,i,j,k,DX);
351
352 //Set::Vector f = ddw(i,j,k)*gradgradu;
353 Set::Vector f = Numeric::Divergence(sig, i, j, k, DX);
354
355 //MATRIX4 AMREX_D_DECL(
356 // Cgrad1 = (Numeric::Stencil<MATRIX4, 1, 0, 0>::D(ddw, i, j, k, 0, DX)),
357 // Cgrad2 = (Numeric::Stencil<MATRIX4, 0, 1, 0>::D(ddw, i, j, k, 0, DX)),
358 // Cgrad3 = (Numeric::Stencil<MATRIX4, 0, 0, 1>::D(ddw, i, j, k, 0, DX)));
359
360 //f += AMREX_D_TERM( ( Cgrad1*gradu).col(0),
361 // +(Cgrad2*gradu).col(1),
362 // +(Cgrad3*gradu).col(2));
363
364 Set::Vector lapv = Numeric::Laplacian(v, i, j, k, DX);
365 f += mu_newton * lapv + b(i,j,k) - mu_dashpot*v(i,j,k);
366
367 Set::Vector udotdot = f / rho;
368 vnew(i, j, k) = v(i, j, k) + dt * udotdot;
369 unew(i, j, k) = u(i, j, k) + dt * v(i, j, k);
370 }
371 });
372 }
373 }
374
375 for (amrex::MFIter mfi(*disp_mf[lev], false); mfi.isValid(); ++mfi)
376 {
377 amrex::Box bx = mfi.grownnodaltilebox() & domain;
378 amrex::Array4<Set::Matrix> const& eps = (*strain_mf[lev]).array(mfi);
379 amrex::Array4<Set::Matrix> const& sig = (*stress_mf[lev]).array(mfi);
380 amrex::Array4<MODEL> const& model = (*model_mf[lev]).array(mfi);
381 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
382 {
383 model(i, j, k).Advance(dt, eps(i, j, k), sig(i, j, k),time);
384 });
385 }
386
387 model_mf[lev]->setMultiGhost(true);
388 model_mf[lev]->FillBoundaryAndSync(geom[lev].periodicity());
389 }
390
391 void Integrate(int amrlev, Set::Scalar /*time*/, int /*step*/,
392 const amrex::MFIter& mfi, const amrex::Box& a_box) override
393 {
394 BL_PROFILE("Integrator::Base::Mechanics::Integrate");
395 if (m_type == Type::Disable) return;
396
397 if (amrex::ParallelDescriptor::NProcs() > 1 && a_box.contains(amrex::IntVect::TheZeroVector()) &&
398 amrlev == 0)
399 {
400 Util::Warning(INFO,"There is a known bug when calculating trac/disp in Base::Mechanics in parallel.");
401 Util::Warning(INFO,"The thermo.dat values likely will not be correct; use the boxlib output instead.");
402 }
403
404 const amrex::Real* DX = geom[amrlev].CellSize();
405 amrex::Box domain = geom[amrlev].Domain();
406 domain.convert(amrex::IntVect::TheNodeVector());
407
408 amrex::Box box = a_box;
409 box.convert(amrex::IntVect::TheNodeVector());
410
411
412 //Set::Scalar dv = AMREX_D_TERM(DX[0], *DX[1], *DX[2]);
413#if AMREX_SPACEDIM == 2
414 Set::Vector da0(DX[1], 0);
415 Set::Vector da1(0, DX[0]);
416#elif AMREX_SPACEDIM == 3
417 Set::Vector da(DX[1] * DX[2], 0, 0);
418#endif
419
420 const Dim3 /*lo= amrex::lbound(domain),*/ hi = amrex::ubound(domain);
421 const Dim3 /*boxlo= amrex::lbound(box),*/ boxhi = amrex::ubound(box);
422
423 amrex::Array4<const Set::Matrix> const& stress = (*stress_mf[amrlev]).array(mfi);
424 amrex::Array4<const Set::Vector> const& disp = (*disp_mf[amrlev]).array(mfi);
425 amrex::ParallelFor(box, [=] AMREX_GPU_DEVICE(int i, int j, int k)
426 {
427#if AMREX_SPACEDIM == 2
428 if (i == hi.x && j < boxhi.y)
429 {
430 trac_hi[0] += (0.5 * (stress(i, j, k) + stress(i, j + 1, k)) * da0);
431 disp_hi[0] = disp(i, j, k);
432 }
433 if (j == hi.y && i < boxhi.x)
434 {
435 trac_hi[1] += (0.5 * (stress(i, j, k) + stress(i + 1, j, k)) * da1);
436 disp_hi[1] = disp(i, j, k);
437 }
438#elif AMREX_SPACEDIM == 3
439 if (i == hi.x && (j < boxhi.y && k < boxhi.z))
440 {
441 trac_hi[0] += (0.25 * (stress(i, j, k) + stress(i, j + 1, k)
442 + stress(i, j, k + 1) + stress(i, j + 1, k + 1)) * da);
443 disp_hi[0] = disp(i, j, k);
444 }
445#endif
446 });
447
448 }
449
450 void TagCellsForRefinement(int lev, amrex::TagBoxArray& a_tags, Set::Scalar /*time*/, int /*ngrow*/) override
451 {
452 BL_PROFILE("Integrator::Base::Mechanics::TagCellsForRefinement");
453 if (m_type == Type::Disable) return;
454
455 Set::Vector DX(geom[lev].CellSize());
456 Set::Scalar DXnorm = DX.lpNorm<2>();
457 for (amrex::MFIter mfi(*strain_mf[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi)
458 {
459 amrex::Box bx = mfi.tilebox();
460 bx.convert(amrex::IntVect::TheCellVector());
461 amrex::Array4<char> const& tags = a_tags.array(mfi);
462 amrex::Array4<Set::Matrix> const& eps = strain_mf[lev]->array(mfi);
463 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
464 {
465 Set::Matrix3 grad = Numeric::NodeGradientOnCell(eps, i, j, k, DX.data());
466 if (grad.norm() * DXnorm > m_elastic_ref_threshold)
467 tags(i, j, k) = amrex::TagBox::SET;
468 });
469 }
470 }
471
472protected:
477 bool psi_on = false;
478
479 int m_interval = 0;
481
487
488 // Only use these if using the "dynamics" option
492 //Set::Field<Set::Matrix4<AMREX_SPACEDIM,MODEL::sym>> ddw_mf;
496
497 //Set::Vector trac_lo[AMREX_SPACEDIM];
498 Set::Vector trac_hi[AMREX_SPACEDIM];
499 Set::Vector disp_hi[AMREX_SPACEDIM];
500
501
503
505
508
510 bool m_print_model = false;
511 bool m_print_residual = false;
512 bool m_time_evolving = false;
514
516
517 bool plot_disp = true;
518 bool plot_stress = true;
519 bool plot_strain = true;
520 bool plot_psi = true;
521 bool plot_rhs = true;
522
523
525
527};
528}
529}
530#endif
#define pp_query_validate(...)
Definition ParmParse.H:103
#define pp_query_default(...)
Definition ParmParse.H:102
#define pp_forbid(...)
Definition ParmParse.H:110
#define pp_queryclass(...)
Definition ParmParse.H:109
#define INFO
Definition Util.H:24
void SetTime(const Set::Scalar a_time)
Definition Elastic.H:49
virtual void Init(amrex::MultiFab *a_rhs, const amrex::Geometry &a_geom, bool a_homogeneous=false) const =0
virtual std::array< Type, AMREX_SPACEDIM > getType(const int &i, const int &j, const int &k, const amrex::Box &domain)=0
Pure abstract IC object from which all other IC objects inherit.
Definition IC.H:23
void Initialize(const int &a_lev, Set::Field< T > &a_field, Set::Scalar a_time=0.0)
Definition IC.H:39
Initialize using a trigonometric series.
Definition Trig.H:23
bool contains(std::string name)
Definition ParmParse.H:173
void select(std::string name, PTRTYPE *&ic_eta, Args &&... args)
Definition ParmParse.H:1056
int query_default(std::string name, T &value, T defaultvalue)
Definition ParmParse.H:293
Set::Field< Set::Vector > disp_mf
Definition Mechanics.H:482
virtual void TimeStepBegin(Set::Scalar a_time, int a_step) override
Definition Mechanics.H:170
Set::Field< Set::Matrix > strain_mf
Definition Mechanics.H:486
Set::Field< MATRIX4 > ddw_mf
Definition Mechanics.H:475
IC::IC< Set::Vector > * ic_rhs
Definition Mechanics.H:502
void Integrate(int amrlev, Set::Scalar, int, const amrex::MFIter &mfi, const amrex::Box &a_box) override
Definition Mechanics.H:391
Set::Field< Set::Vector > vel_old_mf
Definition Mechanics.H:491
Set::Field< Set::Vector > vel_mf
Definition Mechanics.H:490
static void Parse(Mechanics &value, IO::ParmParse &pp)
Definition Mechanics.H:42
Set::Scalar m_elastic_ref_threshold
Definition Mechanics.H:509
IC::IC< Set::Vector > * velocity_ic
Definition Mechanics.H:504
Set::Vector disp_hi[AMREX_SPACEDIM]
Definition Mechanics.H:499
void Initialize(int lev) override
Use the #ic object to initialize::Temp.
Definition Mechanics.H:150
Set::Matrix4< AMREX_SPACEDIM, MODEL::sym > MATRIX4
Definition Mechanics.H:473
void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
Definition Mechanics.H:251
Set::Field< Set::Vector > disp_old_mf
Definition Mechanics.H:489
Solver::Nonlocal::Newton< MODEL > solver
Definition Mechanics.H:506
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar, int) override
Definition Mechanics.H:450
Set::Field< Set::Matrix > stress_mf
Definition Mechanics.H:485
Set::Field< Set::Vector > res_mf
Definition Mechanics.H:484
Set::Field< Set::Vector > rhs_mf
Definition Mechanics.H:483
BC::Operator::Elastic::Elastic * bc
Definition Mechanics.H:507
Set::Field< Set::Scalar > psi_mf
Definition Mechanics.H:476
virtual void UpdateModel(int a_step, Set::Scalar a_time)=0
Set::Field< MODEL > model_mf
Definition Mechanics.H:474
Set::Vector trac_hi[AMREX_SPACEDIM]
Definition Mechanics.H:498
std::vector< amrex::Box > box
Definition Integrator.H:465
amrex::Vector< amrex::Real > dt
Timesteps for each level of refinement.
Definition Integrator.H:394
void RegisterGeneralFab(Set::Field< T > &new_fab, int ncomp, int nghost, bool evolving=true)
Add a templated nodal field.
void RegisterIntegratedVariable(Set::Scalar *integrated_variable, std::string name, bool extensive=true)
Register a variable to be integrated over the spatial domain using the Integrate function.
void SetBC(::BC::Operator::Elastic::Elastic *a_bc)
The different types of Boundary Condtiions are listed in the BC::Operator::Elastic documentation.
Definition Elastic.H:61
void SetUniform(bool a_uniform)
Definition Elastic.H:105
virtual void SetHomogeneous(bool a_homogeneous) override
Definition Elastic.H:48
amrex::Array4< T > Patch(int lev, amrex::MFIter &mfi) const &
Definition Set.H:76
int finest_level
Definition Set.H:67
AMREX_FORCE_INLINE Set::Scalar norm()
Definition Matrix3.H:29
void setPsi(Set::Field< Set::Scalar > &a_psi)
Definition Newton.H:46
void Define(Operator::Elastic< T::sym > &a_op)
Definition Newton.H:30
void compLinearSolverResidual(Set::Field< Set::Vector > &a_res_mf, Set::Field< Set::Vector > &a_u_mf, Set::Field< Set::Vector > &a_b_mf)
Definition Newton.H:502
Set::Scalar solve(const Set::Field< Set::Vector > &a_u_mf, const Set::Field< Set::Vector > &a_b_mf, Set::Field< T > &a_model_mf, Real a_tol_rel, Real a_tol_abs, const char *checkpoint_file=nullptr)
Definition Newton.H:279
Collection of numerical integrator objects.
Definition AllenCahn.H:42
AMREX_FORCE_INLINE Set::Scalar Laplacian(const amrex::Array4< const Set::Scalar > &f, const int &i, const int &j, const int &k, const int &m, const Set::Scalar dx[AMREX_SPACEDIM], std::array< StencilType, AMREX_SPACEDIM > &stencil=DefaultType)
Definition Stencil.H:555
static AMREX_FORCE_INLINE std::array< StencilType, AMREX_SPACEDIM > GetStencil(const int i, const int j, const int k, const amrex::Box domain)
Definition Stencil.H:45
AMREX_FORCE_INLINE Set::Vector NodeGradientOnCell(const amrex::Array4< const Set::Scalar > &f, const int &i, const int &j, const int &k, const Set::Scalar dx[AMREX_SPACEDIM])
Definition Stencil.H:820
AMREX_FORCE_INLINE Set::Vector Gradient(const amrex::Array4< const Set::Scalar > &f, const int &i, const int &j, const int &k, const int &m, const Set::Scalar dx[AMREX_SPACEDIM], std::array< StencilType, AMREX_SPACEDIM > stencil=DefaultType)
Definition Stencil.H:681
AMREX_FORCE_INLINE Set::Vector Divergence(const amrex::Array4< const Set::Matrix > &dw, const int &i, const int &j, const int &k, const Set::Scalar DX[AMREX_SPACEDIM], std::array< StencilType, AMREX_SPACEDIM > &stencil=DefaultType)
Definition Stencil.H:590
AMREX_FORCE_INLINE std::pair< Set::Vector, Set::Matrix > GradientSplit(const amrex::Array4< const Set::Vector > &f, const int &i, const int &j, const int &k, const Set::Scalar dx[AMREX_SPACEDIM], std::array< StencilType, AMREX_SPACEDIM > stencil=DefaultType)
Definition Stencil.H:790
amrex::Real Scalar
Definition Base.H:18
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition Base.H:19
AMREX_FORCE_INLINE Vector Normal(AMREX_D_DECL(bool xmin, bool ymin, bool zmin), AMREX_D_DECL(bool xmax, bool ymax, bool zmax))
Definition Base.H:141
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
Definition Base.H:22
void Abort(const char *msg)
Definition Util.cpp:268
void Warning(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:202