1#ifndef MODEL_SOLID_AFFINE_ISOTROPIC_H_
2#define MODEL_SOLID_AFFINE_ISOTROPIC_H_
24 Define(a_mu, a_lambda, a_F0);
40 return 0.5 * ((
F -
F0).transpose() * (
ddw * ((
F -
F0)))).trace();
50 void Print(std::ostream& out)
const
52 out <<
"lambda=" <<
ddw.Lambda() <<
", mu=" <<
ddw.Mu() <<
", F0=\n" <<
F0;
56 if (
F0.hasNaN())
return true;
57 if (
ddw.contains_nan())
return true;
84 std::pair<std::string, Set::Scalar> moduli[2];
86 pp.
forbid(
"lame",
"Use 'lambda' instead for lame constant");
87 pp.
forbid(
"shear",
"Use 'mu' instead for shear modulus");
88 pp.
forbid(
"bulk",
"Use 'kappa' instead for bulk modulus");
94 pp.
query_exactly<2>({
"lambda",
"mu",
"E",
"nu",
"kappa"}, moduli);
96 if (moduli[0].first ==
"lambda" && moduli[1].first ==
"mu")
98 lambda = moduli[0].second;
99 mu = moduli[1].second;
101 else if (moduli[0].first ==
"lambda" && moduli[1].first ==
"E")
103 lambda = moduli[0].second;
105 mu = (E - 3.0 * lambda + sqrt(E * E + 9.0 * lambda * lambda + 2.0 * E * lambda)) / 4.0;
107 else if (moduli[0].first ==
"lambda" && moduli[1].first ==
"nu")
109 lambda = moduli[0].second;
111 mu = lambda * (1.0 - 2.0 * nu) / 2.0 / nu;
113 else if (moduli[0].first ==
"mu" && moduli[1].first ==
"E")
115 mu = moduli[0].second;
117 lambda = mu * (E - 2.0 * mu) / (3.0 * mu - E);
119 else if (moduli[0].first ==
"mu" && moduli[1].first ==
"nu")
121 mu = moduli[0].second;
123 lambda = 2.0 * mu * nu / (1.0 - 2.0 * nu);
125 else if (moduli[0].first ==
"mu" && moduli[1].first ==
"K")
127 mu = moduli[0].second;
129 lambda = K - (2.0 * mu / 3.0);
131 else if (moduli[0].first ==
"E" && moduli[1].first ==
"nu")
133 Set::Scalar E = moduli[0].second, nu = moduli[1].second;
134 lambda = E * nu / (1.0 + nu) / (1.0 - 2.0 * nu);
135 mu = E / 2.0 / (1.0 + nu);
139 Util::Exception(
INFO,
"Haven't implemented ",moduli[0].first,
" and ",moduli[1].first,
" (sorry!)");
148#define OP_CLASS Isotropic
149#define OP_VARS X(ddw) X(F0)
167 if (
i == 0)
return name +
"_lambda";
168 if (
i == 1)
return name +
"_mu";
170 if (
i == 2)
return name +
"_F0_xx";
171 if (
i == 3)
return name +
"_F0_xy";
172 if (
i == 4)
return name +
"_F0_yx";
173 if (
i == 5)
return name +
"_F0_yy";
174#elif AMREX_SPACEDIM==3
175 if (
i == 2)
return name +
"_F0_xx";
176 if (
i == 3)
return name +
"_F0_xy";
177 if (
i == 4)
return name +
"_F0_xy";
178 if (
i == 5)
return name +
"_F0_yx";
179 if (
i == 6)
return name +
"_F0_yy";
180 if (
i == 7)
return name +
"_F0_yz";
181 if (
i == 8)
return name +
"_F0_zy";
182 if (
i == 9)
return name +
"_F0_zx";
183 if (
i == 10)
return name +
"_F0_zz";
191 for (amrex::MFIter
mfi(
a_dst, amrex::TilingIfNotGPU());
mfi.isValid(); ++
mfi)
193 const amrex::Box&
bx =
mfi.growntilebox(amrex::IntVect(
a_nghost));
196 amrex::Array4<const Model::Solid::Affine::Isotropic>
const&
src = ((*this)[
a_lev])->
array(
mfi);
197 amrex::Array4<Set::Scalar>
const&
dst =
a_dst.array(
mfi);
208#elif AMREX_SPACEDIM==3
#define ALAMO_SINGLE_DEFINITION
bool contains(std::string name)
void forbid(std::string name, std::string explanation)
void query_exactly(std::vector< std::string > names, std::pair< std::string, Set::Scalar > values[N], std::vector< Unit > units=std::vector< Unit >())
Set::Matrix4< AMREX_SPACEDIM, SYM > ddw
Isotropic(Set::Scalar a_mu, Set::Scalar a_lambda, Set::Matrix a_F0=Set::Matrix::Zero())
static void Parse(Isotropic &value, IO::ParmParse &pp)
void Define(Set::Scalar a_mu, Set::Scalar a_lambda, Set::Matrix a_F0)
void Define(Set::Scalar a_mu, Set::Scalar a_lambda)
static Isotropic Random()
void Print(std::ostream &out) const
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::Isotropic > DDW(const Set::Matrix &) const
Set::Scalar W(const Set::Matrix &F) const
Isotropic(Affine< Set::Sym::Isotropic > base)
Set::Matrix DW(const Set::Matrix &F) const
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)