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