Alamo
Isotropic.H
Go to the documentation of this file.
1#ifndef MODEL_SOLID_AFFINE_ISOTROPIC_H_
2#define MODEL_SOLID_AFFINE_ISOTROPIC_H_
3
4#include "AMReX.H"
5#include "IO/ParmParse.H"
7
8namespace Model
9{
10namespace Solid
11{
12namespace Affine
13{
14class Isotropic : public Affine<Set::Sym::Isotropic>
15{
16public:
17
18 Isotropic() = default;
19 ~Isotropic() = default;
20
22 Isotropic(Set::Scalar a_mu, Set::Scalar a_lambda, Set::Matrix a_F0 = Set::Matrix::Zero())
23 {
24 Define(a_mu, a_lambda, a_F0);
25 };
26
27 void Define(Set::Scalar a_mu, Set::Scalar a_lambda, Set::Matrix a_F0)
28 {
29 F0 = a_F0;
31 }
32
33 void Define(Set::Scalar a_mu, Set::Scalar a_lambda)
34 {
36 }
37
38 Set::Scalar W(const Set::Matrix& F) const
39 {
40 return 0.5 * ((F - F0).transpose() * (ddw * ((F - F0)))).trace();
41 }
42 Set::Matrix DW(const Set::Matrix& F) const
43 {
44 return ddw * (F - F0);
45 }
47 {
48 return ddw;
49 }
50 void Print(std::ostream& out) const
51 {
52 out << "lambda=" << ddw.Lambda() << ", mu=" << ddw.Mu() << ", F0=\n" << F0;
53 }
55 {
56 if (F0.hasNaN()) return true;
57 if (ddw.contains_nan()) return true;
58 else return false;
59 }
60public:
61 static Isotropic Zero()
62 {
63 Isotropic ret;
64 Set::Scalar mu = 0.0;
65 Set::Scalar lambda = 0.0;
66 Set::Matrix F0 = Set::Matrix::Zero();
67 ret.Define(mu, lambda, F0);
68 return ret;
69 }
71 {
72 Isotropic ret;
74 Set::Scalar lambda = Util::Random();
75 Set::Matrix F0 = Set::Matrix::Random();
76 ret.Define(mu, lambda, F0);
77 return ret;
78 }
79 static void Parse(Isotropic& value, IO::ParmParse& pp)
80 {
81 Set::Scalar mu = NAN, lambda = NAN;
82 Set::Matrix F0 = Set::Matrix::Zero();
83
84 std::pair<std::string, Set::Scalar> moduli[2];
85
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");
89
90 // Specify exactly two of: lame constant \‍(\lambda\‍),
91 // shear modulus \‍(\mu\‍), Young's modulus \‍(E\‍), Poisson's ratio \‍(\nu\‍),
92 // bulk modulus \‍(\kappa\‍).
93 // \‍(\mu\‍) and \‍(\lambda\‍) are how the final values are stored.
94 pp.query_exactly<2>({"lambda","mu","E","nu","kappa"}, moduli);
95
96 if (moduli[0].first == "lambda" && moduli[1].first == "mu")
97 {
98 lambda = moduli[0].second;
99 mu = moduli[1].second;
100 }
101 else if (moduli[0].first == "lambda" && moduli[1].first == "E")
102 {
103 lambda = moduli[0].second;
104 Set::Scalar E = moduli[1].second;
105 mu = (E - 3.0 * lambda + sqrt(E * E + 9.0 * lambda * lambda + 2.0 * E * lambda)) / 4.0;
106 }
107 else if (moduli[0].first == "lambda" && moduli[1].first == "nu")
108 {
109 lambda = moduli[0].second;
110 Set::Scalar nu = moduli[1].second;
111 mu = lambda * (1.0 - 2.0 * nu) / 2.0 / nu;
112 }
113 else if (moduli[0].first == "mu" && moduli[1].first == "E")
114 {
115 mu = moduli[0].second;
116 Set::Scalar E = moduli[1].second;
117 lambda = mu * (E - 2.0 * mu) / (3.0 * mu - E);
118 }
119 else if (moduli[0].first == "mu" && moduli[1].first == "nu")
120 {
121 mu = moduli[0].second;
122 Set::Scalar nu = moduli[1].second;
123 lambda = 2.0 * mu * nu / (1.0 - 2.0 * nu);
124 }
125 else if (moduli[0].first == "mu" && moduli[1].first == "K")
126 {
127 mu = moduli[0].second;
128 Set::Scalar K = moduli[1].second;
129 lambda = K - (2.0 * mu / 3.0);
130 }
131 else if (moduli[0].first == "E" && moduli[1].first == "nu")
132 {
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);
136 }
137 else
138 {
139 Util::Exception(INFO,"Haven't implemented ",moduli[0].first," and ",moduli[1].first," (sorry!)");
140 }
141
142 if (pp.contains("F0"))
143 {
144 pp_queryarr("F0", F0); // Eigendeformation gradient
145 }
146 value.Define(mu, lambda, F0);
147 }
148#define OP_CLASS Isotropic
149#define OP_VARS X(ddw) X(F0)
151};
153}
154}
155}
156
157template<>
163template<>
166{
167 if (i == 0) return name + "_lambda";
168 if (i == 1) return name + "_mu";
169#if AMREX_SPACEDIM==2
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";
184#endif
185 return name;
186}
187template<>
190{
191 for (amrex::MFIter mfi(a_dst, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
192 {
193 const amrex::Box& bx = mfi.growntilebox(amrex::IntVect(a_nghost));
194 if (bx.ok())
195 {
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);
198 for (int n = 0; n < AMREX_SPACEDIM * AMREX_SPACEDIM; n++)
199 {
200 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
201 dst(i, j, k, a_dstcomp + 0) = src(i, j, k).ddw.Lambda();
202 dst(i, j, k, a_dstcomp + 1) = src(i, j, k).ddw.Mu();
203#if AMREX_SPACEDIM==2
204 dst(i, j, k, a_dstcomp + 2) = src(i, j, k).F0(0, 0);
205 dst(i, j, k, a_dstcomp + 3) = src(i, j, k).F0(0, 1);
206 dst(i, j, k, a_dstcomp + 4) = src(i, j, k).F0(1, 0);
207 dst(i, j, k, a_dstcomp + 5) = src(i, j, k).F0(1, 1);
208#elif AMREX_SPACEDIM==3
209 dst(i, j, k, a_dstcomp + 2) = src(i, j, k).F0(0, 0);
210 dst(i, j, k, a_dstcomp + 3) = src(i, j, k).F0(0, 1);
211 dst(i, j, k, a_dstcomp + 4) = src(i, j, k).F0(0, 2);
212 dst(i, j, k, a_dstcomp + 5) = src(i, j, k).F0(1, 0);
213 dst(i, j, k, a_dstcomp + 6) = src(i, j, k).F0(1, 1);
214 dst(i, j, k, a_dstcomp + 7) = src(i, j, k).F0(1, 2);
215 dst(i, j, k, a_dstcomp + 8) = src(i, j, k).F0(2, 0);
216 dst(i, j, k, a_dstcomp + 9) = src(i, j, k).F0(2, 1);
217 dst(i, j, k, a_dstcomp + 10) = src(i, j, k).F0(2, 2);
218#endif
219 });
220 }
221 }
222 }
223}
224
225
226
227
228#endif
#define pp_queryarr(...)
Definition ParmParse.H:105
#define ALAMO_SINGLE_DEFINITION
Definition Util.H:28
#define INFO
Definition Util.H:23
bool contains(std::string name)
Definition ParmParse.H:173
void forbid(std::string name, std::string explanation)
Definition ParmParse.H:160
void query_exactly(std::vector< std::string > names, std::pair< std::string, Set::Scalar > values[N], std::vector< Unit > units=std::vector< Unit >())
Definition ParmParse.H:1279
Set::Matrix4< AMREX_SPACEDIM, SYM > ddw
Definition Affine.H:37
Isotropic(Set::Scalar a_mu, Set::Scalar a_lambda, Set::Matrix a_F0=Set::Matrix::Zero())
Definition Isotropic.H:22
static void Parse(Isotropic &value, IO::ParmParse &pp)
Definition Isotropic.H:79
void Define(Set::Scalar a_mu, Set::Scalar a_lambda, Set::Matrix a_F0)
Definition Isotropic.H:27
void Define(Set::Scalar a_mu, Set::Scalar a_lambda)
Definition Isotropic.H:33
static Isotropic Random()
Definition Isotropic.H:70
void Print(std::ostream &out) const
Definition Isotropic.H:50
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::Isotropic > DDW(const Set::Matrix &) const
Definition Isotropic.H:46
Set::Scalar W(const Set::Matrix &F) const
Definition Isotropic.H:38
static Isotropic Zero()
Definition Isotropic.H:61
Isotropic(Affine< Set::Sym::Isotropic > base)
Definition Isotropic.H:21
Set::Matrix DW(const Set::Matrix &F) const
Definition Isotropic.H:42
std::string Name(int) const
Definition Set.H:73
void Copy(int, amrex::MultiFab &, int, int) const
Definition Set.H:68
int NComp() const
Definition Set.H:72
A collection of data types and symmetry-reduced data structures.
Definition Base.H:17
amrex::Real Scalar
Definition Base.H:18
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition Base.H:19
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
Definition Base.H:22
Set::Scalar Random()
Definition Set.cpp:34
void Exception(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:225