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 
8 namespace Model
9 {
10 namespace Solid
11 {
12 namespace Affine
13 {
14 class Isotropic : public Affine<Set::Sym::Isotropic>
15 {
16 public:
17 
18  Isotropic() {};
20  Isotropic(Set::Scalar a_mu, Set::Scalar a_lambda, Set::Matrix a_F0 = Set::Matrix::Zero())
21  {
22  Define(a_mu, a_lambda, a_F0);
23  };
24 
25  void Define(Set::Scalar a_mu, Set::Scalar a_lambda, Set::Matrix a_F0)
26  {
27  F0 = a_F0;
29  }
30 
31  void Define(Set::Scalar a_mu, Set::Scalar a_lambda)
32  {
34  }
35 
36  Set::Scalar W(const Set::Matrix& F) const override
37  {
38  return 0.5 * ((F - F0).transpose() * (ddw * ((F - F0)))).trace();
39  }
40  Set::Matrix DW(const Set::Matrix& F) const override
41  {
42  return ddw * (F - F0);
43  }
45  {
46  return ddw;
47  }
48  virtual void Print(std::ostream& out) const override
49  {
50  out << "lambda=" << ddw.Lambda() << ", mu=" << ddw.Mu() << ", F0=\n" << F0;
51  }
52  virtual bool ContainsNan() override
53  {
54  if (F0.hasNaN()) return true;
55  if (ddw.contains_nan()) return true;
56  else return false;
57  }
58 public:
59  static Isotropic Zero()
60  {
61  Isotropic ret;
62  Set::Scalar mu = 0.0;
63  Set::Scalar lambda = 0.0;
64  Set::Matrix F0 = Set::Matrix::Zero();
65  ret.Define(mu, lambda, F0);
66  return ret;
67  }
68  static Isotropic Random()
69  {
70  Isotropic ret;
72  Set::Scalar lambda = Util::Random();
74  ret.Define(mu, lambda, F0);
75  return ret;
76  }
77  static void Parse(Isotropic& value, IO::ParmParse& pp)
78  {
79  Set::Scalar mu = NAN, lambda = NAN;
80  Set::Matrix F0 = Set::Matrix::Zero();
81  if (pp.contains("lame") && pp.contains("shear"))
82  {
83  pp_query("lame", lambda); // Lame modulus
84  pp_query("shear", mu); // Shear modulus
85  }
86  else if (pp.contains("E") && pp.contains("nu"))
87  {
88  Set::Scalar E, nu;
89  pp_query("E", E); // Elastic modulus
90  pp_query("nu", nu); // Poisson's ratio
91  lambda = E * nu / (1.0 + nu) / (1.0 - 2.0 * nu);
92  mu = E / 2.0 / (1.0 + nu);
93  }
94  else if (pp.contains("E") && pp.contains("shear"))
95  {
96  Set::Scalar E;
97  pp_query("E", E); // Young's modulus
98  pp_query("shear", mu); // Shear modulus
99  lambda = mu * (E - 2.0 * mu) / (3.0 * mu - E);
100  }
101  else if (pp.contains("E") && pp.contains("lame"))
102  {
103  Set::Scalar E;
104  pp_query("E", E); // Young's modulus
105  pp_query("lame", lambda); // Lame parameter
106  mu = (E - 3.0 * lambda + sqrt(E * E + 9.0 * lambda * lambda + 2.0 * E * lambda)) / 4.0;
107  }
108  else if (pp.contains("nu") && pp.contains("shear"))
109  {
110  Set::Scalar nu;
111  pp_query("nu", nu); // Poisson's ratio
112  pp_query("shear", mu); // Shear modulus
113  lambda = 2.0 * mu * nu / (1.0 - 2.0 * nu);
114  }
115  else if (pp.contains("lambda") && pp.contains("nu"))
116  {
117  Set::Scalar nu;
118  pp_query("lambda", lambda); // Lame parameter
119  pp_query("nu", nu); // Poisson's ratio
120  mu = lambda * (1.0 - 2.0 * nu) / 2.0 / nu;
121  }
122  else if (pp.contains("bulk") && pp.contains("shear"))
123  {
124  Set::Scalar K;
125  pp_query("bulk", K); // Bulk modulus
126  pp_query("shear", mu); // Shear modulus
127  lambda = K - (2.0 * mu / 3.0);
128  }
129  else
130  {
131  Util::Exception(INFO, "Model parameters not specified with either (lame, shear), or (E, nu)");
132  }
133  if (pp.contains("F0"))
134  {
135  pp_queryarr("F0", F0); // Eigendeformation gradient
136  }
137  value.Define(mu, lambda, F0);
138  }
139 #define OP_CLASS Isotropic
140 #define OP_VARS X(ddw) X(F0)
142 };
144 }
145 }
146 }
147 
148 template<>
151 {
152  return AMREX_SPACEDIM * AMREX_SPACEDIM + 2;
153 }
154 template<>
157 {
158  if (i == 0) return name + "_lambda";
159  if (i == 1) return name + "_mu";
160 #if AMREX_SPACEDIM==2
161  if (i == 2) return name + "_F0_xx";
162  if (i == 3) return name + "_F0_xy";
163  if (i == 4) return name + "_F0_yx";
164  if (i == 5) return name + "_F0_yy";
165 #elif AMREX_SPACEDIM==3
166  if (i == 2) return name + "_F0_xx";
167  if (i == 3) return name + "_F0_xy";
168  if (i == 4) return name + "_F0_xy";
169  if (i == 5) return name + "_F0_yx";
170  if (i == 6) return name + "_F0_yy";
171  if (i == 7) return name + "_F0_yz";
172  if (i == 8) return name + "_F0_zy";
173  if (i == 9) return name + "_F0_zx";
174  if (i == 10) return name + "_F0_zz";
175 #endif
176  return name;
177 }
178 template<>
180 void Set::Field<Model::Solid::Affine::Isotropic>::Copy(int a_lev, amrex::MultiFab& a_dst, int a_dstcomp, int a_nghost) const
181 {
182  for (amrex::MFIter mfi(a_dst, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
183  {
184  const amrex::Box& bx = mfi.growntilebox(amrex::IntVect(a_nghost));
185  if (bx.ok())
186  {
187  amrex::Array4<const Model::Solid::Affine::Isotropic> const& src = ((*this)[a_lev])->array(mfi);
188  amrex::Array4<Set::Scalar> const& dst = a_dst.array(mfi);
189  for (int n = 0; n < AMREX_SPACEDIM * AMREX_SPACEDIM; n++)
190  {
191  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
192  dst(i, j, k, a_dstcomp + 0) = src(i, j, k).ddw.Lambda();
193  dst(i, j, k, a_dstcomp + 1) = src(i, j, k).ddw.Mu();
194 #if AMREX_SPACEDIM==2
195  dst(i, j, k, a_dstcomp + 2) = src(i, j, k).F0(0, 0);
196  dst(i, j, k, a_dstcomp + 3) = src(i, j, k).F0(0, 1);
197  dst(i, j, k, a_dstcomp + 4) = src(i, j, k).F0(1, 0);
198  dst(i, j, k, a_dstcomp + 5) = src(i, j, k).F0(1, 1);
199 #elif AMREX_SPACEDIM==3
200  dst(i, j, k, a_dstcomp + 2) = src(i, j, k).F0(0, 0);
201  dst(i, j, k, a_dstcomp + 3) = src(i, j, k).F0(0, 1);
202  dst(i, j, k, a_dstcomp + 4) = src(i, j, k).F0(0, 2);
203  dst(i, j, k, a_dstcomp + 5) = src(i, j, k).F0(1, 0);
204  dst(i, j, k, a_dstcomp + 6) = src(i, j, k).F0(1, 1);
205  dst(i, j, k, a_dstcomp + 7) = src(i, j, k).F0(1, 2);
206  dst(i, j, k, a_dstcomp + 8) = src(i, j, k).F0(2, 0);
207  dst(i, j, k, a_dstcomp + 9) = src(i, j, k).F0(2, 1);
208  dst(i, j, k, a_dstcomp + 10) = src(i, j, k).F0(2, 2);
209 #endif
210  });
211  }
212  }
213  }
214 }
215 
216 
217 
218 
219 #endif
Set::Sym
Sym
Definition: Base.H:197
Model::Solid::Affine::Isotropic::DDW
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::Isotropic > DDW(const Set::Matrix &) const override
Definition: Isotropic.H:44
Model::Solid::Affine::Isotropic::Define
void Define(Set::Scalar a_mu, Set::Scalar a_lambda, Set::Matrix a_F0)
Definition: Isotropic.H:25
Model::Solid::Affine::Isotropic::Isotropic
Isotropic(Affine< Set::Sym::Isotropic > base)
Definition: Isotropic.H:19
Model::Solid::Affine::Isotropic
Definition: Isotropic.H:14
Model::Solid::Affine::Isotropic::ContainsNan
virtual bool ContainsNan() override
Definition: Isotropic.H:52
ALAMO_SINGLE_DEFINITION
#define ALAMO_SINGLE_DEFINITION
Definition: Util.H:25
Model::Solid::Affine::Isotropic::Zero
static Isotropic Zero()
Definition: Isotropic.H:59
Set::Field::Copy
void Copy(int, amrex::MultiFab &, int, int) const
Definition: Set.H:68
Affine.H
ParmParse.H
Model::Solid::Affine::Affine< Set::Sym::Isotropic >::F0
Set::Matrix F0
Definition: Affine.H:38
pp_query
#define pp_query(...)
Definition: ParmParse.H:104
Model::Solid::Affine::Isotropic::W
Set::Scalar W(const Set::Matrix &F) const override
Definition: Isotropic.H:36
Model::Solid::Affine::Isotropic::Print
virtual void Print(std::ostream &out) const override
Definition: Isotropic.H:48
Set::Field::Name
std::string Name(int) const
Definition: Set.H:73
Util::Random
Set::Scalar Random()
Definition: Set.cpp:9
Set::Scalar
amrex::Real Scalar
Definition: Base.H:19
Set
A collection of data types and symmetry-reduced data structures.
Definition: Base.H:17
InClassOperators.H
Set::Matrix
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
Definition: Base.H:23
ExtClassOperators.H
Model::Solid::Affine::Isotropic::Parse
static void Parse(Isotropic &value, IO::ParmParse &pp)
Definition: Isotropic.H:77
Model::Solid::Affine::Isotropic::Define
void Define(Set::Scalar a_mu, Set::Scalar a_lambda)
Definition: Isotropic.H:31
Util::Exception
void Exception(std::string file, std::string func, int line, Args const &... args)
Definition: Util.H:194
Model::Solid::Affine::Isotropic::Isotropic
Isotropic(Set::Scalar a_mu, Set::Scalar a_lambda, Set::Matrix a_F0=Set::Matrix::Zero())
Definition: Isotropic.H:20
IO::ParmParse::contains
bool contains(std::string name)
Definition: ParmParse.H:151
Model::Solid::Affine::Isotropic::Random
static Isotropic Random()
Definition: Isotropic.H:68
Model::Solid::Affine::Isotropic::Isotropic
Isotropic()
Definition: Isotropic.H:18
Model::Solid::Affine::Isotropic::DW
Set::Matrix DW(const Set::Matrix &F) const override
Definition: Isotropic.H:40
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::Isotropic >
Model::Solid::Affine::Affine
Definition: Affine.H:34
Model::Solid::Affine::Affine< Set::Sym::Isotropic >::ddw
Set::Matrix4< AMREX_SPACEDIM, SYM > ddw
Definition: Affine.H:37
IO::ParmParse
Definition: ParmParse.H:110
Set::Field::NComp
int NComp() const
Definition: Set.H:72
INFO
#define INFO
Definition: Util.H:20
Model
Definition: Constant.H:12
Model::Solid::F
@ F
Definition: Solid.H:26
pp_queryarr
#define pp_queryarr(...)
Definition: ParmParse.H:103