Alamo
Mechanics.H
Go to the documentation of this file.
1//
2// This is a general purpose integrator that focuses on solving
3// elasticity/mechanics equations in the absence of other multiphysics
4// simulations.
5// It is enabled by :code:`alamo.program=mechanics` if used on its own, in
6// which case there is no prefix.
7// If it is being used by another integrator, see that integrator to determine
8// the value of :code:`[prefix]` (often equal to :code:`elastic`).
9//
10// This integrator inherits from :ref:`Integrator::Base::Mechanics`; see
11// documentation for that integrator for additional parameters.
12//
13// :bdg-primary-line:`Model setup`
14// There are two basic tools for setting up a mechanics problem in this
15// integrator.
16//
17// 1. The :code:`eta` field: this is used to mix models of different types.
18// Use :code:`nmodels` to specify how many material models to use, and then
19// specify each model as :code:`model1`, :code:`model2`, etc.
20// The type of moel is specified using the :code:`alamo.program.mechanics.model`
21// input.
22//
23// Once you have done this, you must determine the spatial distribution of each
24// model. This is done by choosing an IC for the eta field with :code:`ic.type`.
25// The :code:`IC::Expression` is the most general and recommended.
26// The models are then mixed linearly, i.e.
27//
28// .. math::
29//
30// W_{\textrm{eff}} = \sum_{i=1}^N W_i\,\eta_i(\mathbf{x})
31//
32// See the :ref:`Eshelby` test for an example of model mixing.
33//
34// 2. The :code:`psi` field: this is used specifically for cases where a "void" region
35// is desired. Its usage is similar to the :code:`eta` case, and is conceptually
36// similar in that it scales the model field to near-zero in order to mimic the
37// (lack of) mechanical behavior in an empty region.
38// It is important to use :code:`psi` here, for reasons that are discussed in detail
39// in
40// `this paper <https://doi.org/10.1007/s00466-023-02325-8>`_.
41// The initialization of :code:`psi` is similar to that for :code:`eta`.
42//
43// See the :ref:`PlateHole` and :ref:`RubberPlateHole` for canonical exmaples.
44// The :ref:`Integrator::Fracture` and :ref:`Integrator::TopOp` integrators are examples
45// of integrators that leverage the psi property.
46//
47// :bdg-primary-line:`Body forces`
48// currently have limited support due to the relatively low number of
49// times they are needed. See the :ref:`Integrator::Base::Mechanics` documentation for
50// detail.
51// See the :ref:`TrigTest` test for examples of current body force implementation.
52//
53// :bdg-primary-line:`Boundary conditions`
54// are implemented using the :ref:`BC::Operator::Elastic` classes.
55// See the documentation on these classes for more detail.
56// See any of the mechanics-based tests for examples of boundary condition application.
57//
58
59#ifndef INTEGRATOR_MECHANICS_H
60#define INTEGRATOR_MECHANICS_H
61#include <iostream>
62#include <fstream>
63#include <iomanip>
64#include <numeric>
65
66#include "AMReX.H"
67#include "AMReX_ParallelDescriptor.H"
68#include "AMReX_ParmParse.H"
69
70#include "IO/ParmParse.H"
72
73
74#include "IC/IC.H"
75#include "BC/BC.H"
79
80#include "IC/Ellipse.H"
81#include "IC/Voronoi.H"
82#include "IC/Constant.H"
83#include "IC/Expression.H"
84#include "IC/BMP.H"
85#include "IC/PNG.H"
86#include "IC/PSRead.H"
87#include "IC/Constant.H"
88#include "BC/Constant.H"
89#include "Numeric/Stencil.H"
90
91#include "Model/Solid/Solid.H"
94
95#include "Operator/Operator.H"
96#include "Base/Mechanics.H"
97
98namespace Integrator
99{
100template<class MODEL>
101class Mechanics : virtual public Base::Mechanics<MODEL>
102{
103public:
104 static constexpr const char *name = "mechanics";
105
106 Mechanics() : Base::Mechanics<MODEL>() { }
107 Mechanics(IO::ParmParse& pp) : Base::Mechanics<MODEL>()
108 {
109 Parse(*this, pp);
110 }
111
113 {
114 delete ic_eta;
115 delete ic_psi;
116 delete ic_trac_normal;
117 delete bc_psi;
118 delete bc_trac_normal;
119 }
120
121 // Mechanics inputs. See also :ref:`Integrator::Base::Mechanics`
122 static void Parse(Mechanics& value, IO::ParmParse& pp)
123 {
125
126 int nmodels;
127 pp_query_default("nmodels", nmodels, 1); // Number of elastic model varieties
128
129 pp.queryclass_enumerate<MODEL>("model",value.models,nmodels);
130
131 Util::Assert(INFO, TEST((int)nmodels == (int)value.models.size()),
132 "nmodels does not match the # of specified models");
133
134 Util::Assert(INFO, TEST(value.models.size() > 0));
135 value.RegisterNodalFab(value.eta_mf, value.models.size(), 2, "eta", true);
136 // Refinement threshold for eta field
137 pp_query_default("eta_ref_threshold", value.m_eta_ref_threshold, 0.01);
138 // Refinement threshold for strain gradient
139 pp_query_default("ref_threshold", value.m_elastic_ref_threshold, 0.01);
140 // Explicity impose neumann condition on model at domain boundaries (2d only)
141 pp_query_default("model_neumann_boundary",value.model_neumann_boundary,false);
142
143 // Read in IC for eta
144 if (nmodels > 1 && pp.contains("ic.type"))
145 {
146 // Select the initial condition for eta
148
149 value.eta_reset_on_regrid = true;
150 // Whether to re-initialize eta when re-gridding occurs.
151 // Default is false unless eta ic is set, then default is.
152 // true.
153 pp_query_default("eta.reset_on_regrid", value.eta_reset_on_regrid, true);
154 }
155
156 // Read in IC for psi
157 if (pp.contains("psi.ic.type"))
158 {
159 // Select initial condition for psi field
161
162 value.bc_psi = new BC::Nothing();
163 value.RegisterNewFab(value.psi_mf, value.bc_psi, 1, 2, "psi", value.plot_psi);
164 value.psi_on = true;
165
166 value.psi_reset_on_regrid = true;
167 // Whether to re-initialize psi when re-gridding occurs.
168 // Default is false unless a psi ic is set, then default is
169 // true.
170 pp_query_default("psi.reset_on_regrid", value.psi_reset_on_regrid, true);
171 }
172
173 if (pp.contains("trac_normal.ic.type"))
174 {
175 // Read in IC for the "normal traction" field (applied at diffuse boundaries)
176 pp.select<IC::Ellipse,IC::Constant,IC::Expression,IC::PSRead>("trac_normal.ic",value.ic_trac_normal,value.geom);
177
178 if (value.ic_trac_normal)
179 {
180 value.bc_trac_normal = new BC::Nothing();
181 value.RegisterNewFab(value.trac_normal_mf, value.bc_trac_normal, 1, 2, "trac_normal", true);
182 }
183 }
185
186 }
187
188 void Initialize(int lev) override
189 {
191 eta_mf[lev]->setVal(0.0);
192 if (models.size() > 1 && ic_eta) ic_eta->Initialize(lev, eta_mf);
193 else eta_mf[lev]->setVal(1.0);
194
195 if (psi_on) ic_psi->Initialize(lev, psi_mf);
197 }
198
199 virtual void UpdateModel(int a_step, Set::Scalar time) override
200 {
202
203 if (ic_trac_normal)
204 {
205 for (int lev = 0; lev <= finest_level; ++lev)
206 {
208 psi_mf[lev]->FillBoundary();
209 amrex::Box domain = this->geom[lev].Domain();
210 domain.convert(amrex::IntVect::TheNodeVector());
211 Set::Vector DX(geom[lev].CellSize());
212
213 for (MFIter mfi(*model_mf[lev], false); mfi.isValid(); ++mfi)
214 {
215 amrex::Box bx = mfi.nodaltilebox();
216 bx = bx & domain;
217
218 amrex::Array4<Set::Vector> const& rhs = rhs_mf[lev]->array(mfi);
219 amrex::Array4<const Set::Scalar> const& psi = psi_mf[lev]->array(mfi);
220 amrex::Array4<const Set::Scalar> const& trac_normal = trac_normal_mf[lev]->array(mfi);
221
222 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
223 Set::Vector grad = Numeric::CellGradientOnNode(psi, i, j, k, 0, DX.data());
224 rhs(i,j,k) = trac_normal(i,j,k) * grad;
225 });
226 }
227 Util::RealFillBoundary(*rhs_mf[lev], geom[lev]);
228 }
229 }
230
231 if (a_step > 0) return;
232
233 for (int lev = 0; lev <= finest_level; ++lev)
234 {
235 eta_mf[lev]->FillBoundary();
236
237 amrex::Box domain = this->geom[lev].Domain();
238 domain.convert(amrex::IntVect::TheNodeVector());
239
240 Set::Vector DX(geom[lev].CellSize());
241
242 for (MFIter mfi(*model_mf[lev], false); mfi.isValid(); ++mfi)
243 {
244 amrex::Box bx = mfi.grownnodaltilebox();
245
246 amrex::Array4<MODEL> const& model = model_mf[lev]->array(mfi);
247 amrex::Array4<const Set::Scalar> const& eta = eta_mf[lev]->array(mfi);
248
249 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
250 model(i, j, k) = MODEL::Zero();
251 for (unsigned int n = 0; n < models.size(); n++)
252 model(i, j, k) += eta(i, j, k, n) * models[n];
253 });
254 }
255
256
258 {
259 Util::AssertException(INFO,TEST(AMREX_SPACEDIM==2),"neumann boundary works in 2d only");
260 for (MFIter mfi(*model_mf[lev], false); mfi.isValid(); ++mfi)
261 {
262 amrex::Box bx = mfi.grownnodaltilebox() & domain;
263 amrex::Array4<MODEL> const& model = this->model_mf[lev]->array(mfi);
264 const Dim3 lo = amrex::lbound(domain), hi = amrex::ubound(domain);
265
266 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
267 {
268 if (i==lo.x && j==lo.y)
269 model(i,j,k) = 0.5*(model(i+1,j,k)+model(i,j+1,k));
270 else if (i==lo.x && j==hi.y)
271 model(i,j,k) = 0.5*(model(i+1,j,k)+model(i,j-1,k));
272 else if (i==hi.x && j==lo.y)
273 model(i,j,k) = 0.5*(model(i-1,j,k)+model(i,j+1,k));
274 else if (i==hi.x && j==hi.y)
275 model(i,j,k) = 0.5*(model(i-1,j,k)+model(i,j-1,k));
276
277 else if (i==lo.x)
278 model(i,j,k) = model(i+1,j,k);
279 else if (i==hi.x)
280 model(i,j,k) = model(i-1,j,k);
281 else if (j==lo.y)
282 model(i,j,k) = model(i,j+1,k);
283 else if (j==hi.y)
284 model(i,j,k) = model(i,j-1,k);
285 });
286 }
287 }
288
289 Util::RealFillBoundary(*model_mf[lev], geom[lev]);
290 }
291 }
292
293 void TagCellsForRefinement(int lev, amrex::TagBoxArray& a_tags, Set::Scalar a_time, int a_ngrow) override
294 {
296 Base::Mechanics<MODEL>::TagCellsForRefinement(lev, a_tags, a_time, a_ngrow);
297
298 Set::Vector DX(geom[lev].CellSize());
299 Set::Scalar DXnorm = DX.lpNorm<2>();
300 for (amrex::MFIter mfi(*eta_mf[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi)
301 {
302 amrex::Box bx = mfi.nodaltilebox();
303 amrex::Array4<char> const& tags = a_tags.array(mfi);
304 amrex::Array4<Set::Scalar> const& eta = eta_mf[lev]->array(mfi);
305 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
306 {
307 auto sten = Numeric::GetStencil(i, j, k, bx);
308 for (int n = 0; n < eta_mf[lev]->nComp(); n++)
309 {
310 Set::Vector grad = Numeric::Gradient(eta, i, j, k, n, DX.data(), sten);
311 if (grad.lpNorm<2>() * DXnorm > m_eta_ref_threshold)
312 tags(i, j, k) = amrex::TagBox::SET;
313 }
314 });
315 if (psi_on)
316 {
317 amrex::Array4<Set::Scalar> const& psi = psi_mf[lev]->array(mfi);
318 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
319 {
320 auto sten = Numeric::GetStencil(i, j, k, bx);
321 {
322 Set::Vector gradpsi = Numeric::Gradient(psi, i, j, k, 0, DX.data(), sten);
323 if (gradpsi.lpNorm<2>() * DXnorm > m_eta_ref_threshold)
324 tags(i, j, k) = amrex::TagBox::SET;
325 }
326 });
327 }
328 }
329 }
330
331 void Regrid(int lev, Set::Scalar time) override
332 {
333 if (eta_reset_on_regrid && models.size() > 1 && ic_eta) ic_eta->Initialize(lev, eta_mf, time);
335 }
336
337protected:
340 std::vector<MODEL> models;
348
350
352
353 using Base::Mechanics<MODEL>::m_type;
354 using Base::Mechanics<MODEL>::finest_level;
355 using Base::Mechanics<MODEL>::geom;
356 using Base::Mechanics<MODEL>::model_mf;
357 using Base::Mechanics<MODEL>::psi_mf;
358 using Base::Mechanics<MODEL>::psi_on;
359 using Base::Mechanics<MODEL>::rhs_mf;
360};
361
362
363
364
365
366
367
368
369
370
371} // namespace Integrator
372#endif
#define pp_query_default(...)
Definition ParmParse.H:102
#define TEST(x)
Definition Util.H:25
#define INFO
Definition Util.H:24
Definition BC.H:43
Definition BMP.H:22
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
Definition PNG.H:26
bool contains(std::string name)
Definition ParmParse.H:173
void select(std::string name, PTRTYPE *&ic_eta, Args &&... args)
Definition ParmParse.H:1056
void queryclass(std::string name, T *value)
Definition ParmParse.H:960
int queryclass_enumerate(std::string a_name, std::vector< T > &value, int number=1)
Definition ParmParse.H:764
Set::Scalar m_elastic_ref_threshold
Definition Mechanics.H:509
void Initialize(int lev) override
Use the #ic object to initialize::Temp.
Definition Mechanics.H:150
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar, int) override
Definition Mechanics.H:450
Set::Field< Set::Vector > rhs_mf
Definition Mechanics.H:483
Set::Field< Set::Scalar > psi_mf
Definition Mechanics.H:476
Set::Field< MODEL > model_mf
Definition Mechanics.H:474
void RegisterNodalFab(Set::Field< Set::Scalar > &new_fab, int ncomp, int nghost, std::string name, bool writeout, bool evolving=true, std::vector< std::string > suffix={})
Add a new node-based scalar field.
void RegisterNewFab(Set::Field< Set::Scalar > &new_fab, BC::BC< Set::Scalar > *new_bc, int ncomp, int nghost, std::string name, bool writeout, bool evolving=true, std::vector< std::string > suffix={})
Add a new cell-based scalar field.
IC::IC< Set::Scalar > * ic_trac_normal
Definition Mechanics.H:343
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar a_time, int a_ngrow) override
Definition Mechanics.H:293
IC::IC< Set::Scalar > * ic_eta
Definition Mechanics.H:341
Set::Field< Set::Scalar > eta_mf
Definition Mechanics.H:338
static void Parse(Mechanics &value, IO::ParmParse &pp)
Definition Mechanics.H:122
void Regrid(int lev, Set::Scalar time) override
Definition Mechanics.H:331
Mechanics(IO::ParmParse &pp)
Definition Mechanics.H:107
std::vector< MODEL > models
Definition Mechanics.H:340
BC::BC< Set::Scalar > * bc_trac_normal
Definition Mechanics.H:345
void Initialize(int lev) override
Definition Mechanics.H:188
BC::BC< Set::Scalar > * bc_psi
Definition Mechanics.H:344
Set::Scalar m_eta_ref_threshold
Definition Mechanics.H:339
static constexpr const char * name
Definition Mechanics.H:104
IC::IC< Set::Scalar > * ic_psi
Definition Mechanics.H:342
Set::Field< Set::Scalar > trac_normal_mf
Definition Mechanics.H:351
virtual void UpdateModel(int a_step, Set::Scalar time) override
Definition Mechanics.H:199
Collection of numerical integrator objects.
Definition AllenCahn.H:42
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 CellGradientOnNode(const amrex::Array4< const Set::Scalar > &f, const int &i, const int &j, const int &k, const int &m, const Set::Scalar dx[AMREX_SPACEDIM])
Definition Stencil.H:699
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::Real Scalar
Definition Base.H:18
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition Base.H:19
AMREX_FORCE_INLINE void AssertException(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)
Definition Util.H:254
AMREX_FORCE_INLINE void Assert(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)
Definition Util.H:58
void Message(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:129
AMREX_FORCE_INLINE void RealFillBoundary(amrex::FabArray< amrex::BaseFab< T > > &a_mf, const amrex::Geometry &, const int nghost=2)
Definition Util.H:339