52#ifndef MODEL_SOLID_LINEAR_ISOTROPIC_H_
53#define MODEL_SOLID_LINEAR_ISOTROPIC_H_
94 void Print(std::ostream &out)
const
119 bool planestress =
false;
120 std::pair<std::string, Set::Scalar> moduli[2];
122 pp.
forbid(
"lame",
"Use 'lambda' instead for lame constant");
123 pp.
forbid(
"shear",
"Use 'mu' instead for shear modulus");
124 pp.
forbid(
"bulk",
"Use 'kappa' instead for bulk modulus");
130 pp.
query_exactly<2>({
"lambda",
"mu",
"E",
"nu",
"kappa"}, moduli);
132 if (moduli[0].first ==
"lambda" && moduli[1].first ==
"mu")
134 lambda = moduli[0].second;
135 mu = moduli[1].second;
137 else if (moduli[0].first ==
"lambda" && moduli[1].first ==
"E")
139 lambda = moduli[0].second;
141 mu = (E - 3.0 * lambda + sqrt(E * E + 9.0 * lambda * lambda + 2.0 * E * lambda)) / 4.0;
143 else if (moduli[0].first ==
"lambda" && moduli[1].first ==
"nu")
145 lambda = moduli[0].second;
147 mu = lambda * (1.0 - 2.0 * nu) / 2.0 / nu;
149 else if (moduli[0].first ==
"mu" && moduli[1].first ==
"E")
151 mu = moduli[0].second;
153 lambda = mu * (E - 2.0 * mu) / (3.0 * mu - E);
155 else if (moduli[0].first ==
"mu" && moduli[1].first ==
"nu")
157 mu = moduli[0].second;
159 lambda = 2.0 * mu * nu / (1.0 - 2.0 * nu);
161 else if (moduli[0].first ==
"mu" && moduli[1].first ==
"K")
163 mu = moduli[0].second;
165 lambda = K - (2.0 * mu / 3.0);
167 else if (moduli[0].first ==
"E" && moduli[1].first ==
"nu")
169 Set::Scalar E = moduli[0].second, nu = moduli[1].second;
170 lambda = E * nu / (1.0 + nu) / (1.0 - 2.0 * nu);
171 mu = E / 2.0 / (1.0 + nu);
175 Util::Exception(
INFO,
"Haven't implemented ",moduli[0].first,
" and ",moduli[1].first,
" (sorry!)");
183 if (AMREX_SPACEDIM==2 && planestress)
184 value.
Define(mu,lambda*(1.0 - lambda/(2.*mu + lambda)));
189 #define OP_CLASS Isotropic
190 #define OP_VARS X(ddw)
211 if (
i == 0)
return name +
"_lambda";
212 if (
i == 1)
return name +
"_mu";
219 for (amrex::MFIter
mfi(
a_dst, amrex::TilingIfNotGPU());
mfi.isValid(); ++
mfi)
221 const amrex::Box&
bx =
mfi.growntilebox(amrex::IntVect(
a_nghost));
224 amrex::Array4<const Model::Solid::Linear::Isotropic>
const&
src = ((*this)[
a_lev])->
array(
mfi);
225 amrex::Array4<Set::Scalar>
const&
dst =
a_dst.array(
mfi);
#define ALAMO_SINGLE_DEFINITION
void forbid(std::string name, std::string explanation)
int query_default(std::string name, T &value, T defaultvalue)
void query_exactly(std::vector< std::string > names, std::pair< std::string, Set::Scalar > values[N], std::vector< Unit > units=std::vector< Unit >())
Set::Matrix DW(const Set::Matrix &gradu) const
Isotropic(Solid< Set::Sym::Isotropic > base)
Isotropic(Set::Scalar a_mu, Set::Scalar a_lambda)
static const KinematicVariable kinvar
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::Isotropic > DDW(const Set::Matrix &) const
void Print(std::ostream &out) const
Set::Scalar W(const Set::Matrix &gradu) const
static Isotropic Random()
static void Parse(Isotropic &value, IO::ParmParse &pp)
void Define(Set::Scalar a_mu, Set::Scalar a_lambda)
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::Isotropic > ddw
std::string Name(int) const
void Copy(int, amrex::MultiFab &, int, int) const
A collection of data types and symmetry-reduced data structures.
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
void Exception(std::string file, std::string func, int line, Args const &... args)