Line data Source code
1 : //
2 : // This model implements an isotropic linear elastic material.
3 : // See `this link <https://en.wikipedia.org/wiki/Linear_elasticity#(An)isotropic_(in)homogeneous_media>`_
4 : // for more information about the theory.
5 : //
6 : // Free energy for a linear material is defined as
7 : //
8 : // .. math::
9 : //
10 : // W(\nabla\mathbf{u}) =
11 : // \frac{1}{2}\nabla\mathbf{u}\cdot\mathbb{C}\,\nabla\mathbf{u}
12 : //
13 : // For an isotropic material, stress and strain are related through
14 : //
15 : // .. math::
16 : //
17 : // \mathbb{C}_{ijkl} = \lambda \delta_{ij}\varepsilon_{kk} + 2\mu\varepsilon_{ij}
18 : //
19 : // where :math:`\lambda` and :math:`\mu` are the Lame constant and shear modulus, respectively.
20 : // Users can specify these either through (:code:`lame` and :code:`shear`)
21 : // OR (:code:`lambda` and :code:`mu`) OR (:code:`E` and :code:`nu`).
22 : //
23 : // Class methods:
24 : //
25 : // #. :code:`Isotropic()`:
26 : // Basic constructor. Does nothing, and leaves all values initiated as NAN.
27 : // #. :code:`Isotropic(Solid<Set::Sym::Isotropic> base)`
28 : // Basic constructor. Does nothing gut allows for inheritance.
29 : // #. :code:`Isotropic(Set::Scalar a_mu, Set::Scalar a_lambda)`
30 : // BAD old-fashioned constructor. Do not use!
31 : // #. :code:`~Isotropic()`
32 : // Simple destructor. Don't need to change it.
33 : // #. :code:`void Define(Set::Scalar a_mu, Set::Scalar a_lambda)`
34 : // BAD old-fashioned way of doing things. Use :code:`Parse` instead.
35 : // #. :code:`Set::Scalar W(const Set::Matrix & gradu) const override`
36 : // Returns elastic free energy density
37 : // #. :code:`Set::Matrix DW(const Set::Matrix & gradu) const override`
38 : // Returns first derivative of free energy, the stress tensor
39 : // #. :code:`Set::Matrix4<...> DDW(const Set::Matrix & ) const override`
40 : // Returns second derivative of free energy, the modulus tensor
41 : // #. :code:`virtual void Print(std::ostream &out) const override`
42 : // Prints the modulus tensor object to output stream (usually the terminal)
43 : // #. :code:`static Isotropic Random()`
44 : // Static method that generates a random yet acceptable model.
45 : // #. :code:`static Isotropic Zero()`
46 : // Static method that generates a "zero" element (so that there is no effect under addition)
47 : // #. :code:`static void Parse(Isotropic & value, IO::ParmParse & pp)`
48 : // Parser where all the IO occurs
49 : //
50 : //
51 :
52 : #ifndef MODEL_SOLID_LINEAR_ISOTROPIC_H_
53 : #define MODEL_SOLID_LINEAR_ISOTROPIC_H_
54 :
55 : #include "Model/Solid/Solid.H"
56 : #include "IO/ParmParse.H"
57 :
58 : namespace Model
59 : {
60 : namespace Solid
61 : {
62 : namespace Linear
63 : {
64 : class Isotropic : public Solid<Set::Sym::Isotropic>
65 : {
66 : public:
67 :
68 207894 : Isotropic() {};
69 : Isotropic(Solid<Set::Sym::Isotropic> base) : Solid<Set::Sym::Isotropic>(base) {};
70 : Isotropic(Set::Scalar a_mu, Set::Scalar a_lambda)
71 : {
72 : Define(a_mu,a_lambda);
73 : };
74 283261 : virtual ~Isotropic() {};
75 :
76 49906 : void Define(Set::Scalar a_mu, Set::Scalar a_lambda)
77 : {
78 49906 : ddw = Set::Matrix4<AMREX_SPACEDIM,Set::Sym::Isotropic>(a_lambda,a_mu);
79 49906 : }
80 :
81 260 : Set::Scalar W(const Set::Matrix & gradu) const override
82 : {
83 520 : return ( 0.5 * gradu.transpose() * (ddw*gradu) ).trace();
84 : }
85 134932 : Set::Matrix DW(const Set::Matrix & gradu) const override
86 : {
87 269864 : return ddw*gradu;
88 : }
89 66506 : Set::Matrix4<AMREX_SPACEDIM,Set::Sym::Isotropic> DDW(const Set::Matrix & /*gradu*/) const override
90 : {
91 66506 : return ddw;
92 : }
93 0 : virtual void Print(std::ostream &out) const override
94 : {
95 0 : out << ddw;
96 0 : }
97 :
98 : public:
99 : Set::Matrix4<AMREX_SPACEDIM,Set::Sym::Isotropic> ddw;
100 : static const KinematicVariable kinvar = KinematicVariable::gradu;
101 :
102 : public:
103 44 : static Isotropic Random()
104 : {
105 44 : Isotropic ret;
106 44 : ret.Define(Util::Random(),Util::Random());
107 44 : return ret;
108 0 : }
109 49855 : static Isotropic Zero()
110 : {
111 49855 : Isotropic ret;
112 49855 : ret.Define(0.,0.);
113 49855 : return ret;
114 : }
115 7 : static void Parse(Isotropic & value, IO::ParmParse & pp)
116 : {
117 7 : Set::Scalar mu = NAN, lambda = NAN;
118 7 : bool planestress = false;
119 42 : std::pair<std::string, Set::Scalar> moduli[2];
120 :
121 56 : pp.forbid("lame","Use 'lambda' instead for lame constant");
122 56 : pp.forbid("shear","Use 'mu' instead for shear modulus");
123 49 : pp.forbid("bulk","Use 'kappa' instead for bulk modulus");
124 :
125 : // Specify exactly two of: lame constant \(\lambda\),
126 : // shear modulus \(\mu\), Young's modulus \(E\), Poisson's ratio \(\nu\),
127 : // bulk modulus \(\kappa\).
128 : // \(\mu\) and \(\lambda\) are how the final values are stored.
129 14 : pp.query_exactly<2>({"lambda","mu","E","nu","kappa"}, moduli);
130 :
131 7 : if (moduli[0].first == "lambda" && moduli[1].first == "mu")
132 : {
133 0 : lambda = moduli[0].second;
134 0 : mu = moduli[1].second;
135 : }
136 7 : else if (moduli[0].first == "lambda" && moduli[1].first == "E")
137 : {
138 0 : lambda = moduli[0].second;
139 0 : Set::Scalar E = moduli[1].second;
140 0 : mu = (E - 3.0 * lambda + sqrt(E * E + 9.0 * lambda * lambda + 2.0 * E * lambda)) / 4.0;
141 : }
142 7 : else if (moduli[0].first == "lambda" && moduli[1].first == "nu")
143 : {
144 0 : lambda = moduli[0].second;
145 0 : Set::Scalar nu = moduli[1].second;
146 0 : mu = lambda * (1.0 - 2.0 * nu) / 2.0 / nu;
147 : }
148 7 : else if (moduli[0].first == "mu" && moduli[1].first == "E")
149 : {
150 0 : mu = moduli[0].second;
151 0 : Set::Scalar E = moduli[1].second;
152 0 : lambda = mu * (E - 2.0 * mu) / (3.0 * mu - E);
153 : }
154 7 : else if (moduli[0].first == "mu" && moduli[1].first == "nu")
155 : {
156 0 : mu = moduli[0].second;
157 0 : Set::Scalar nu = moduli[1].second;
158 0 : lambda = 2.0 * mu * nu / (1.0 - 2.0 * nu);
159 : }
160 7 : else if (moduli[0].first == "mu" && moduli[1].first == "K")
161 : {
162 0 : mu = moduli[0].second;
163 0 : Set::Scalar K = moduli[1].second;
164 0 : lambda = K - (2.0 * mu / 3.0);
165 : }
166 7 : else if (moduli[0].first == "E" && moduli[1].first == "nu")
167 : {
168 7 : Set::Scalar E = moduli[0].second, nu = moduli[1].second;
169 7 : lambda = E * nu / (1.0 + nu) / (1.0 - 2.0 * nu);
170 7 : mu = E / 2.0 / (1.0 + nu);
171 : }
172 : else
173 : {
174 0 : Util::Exception(INFO,"Haven't implemented ",moduli[0].first," and ",moduli[1].first," (sorry!)");
175 : }
176 :
177 : // Whether or not to use the
178 : // `plane stress <https://en.wikipedia.org/wiki/Plane_stress>`_
179 : // approximation.
180 36 : pp.query_default("planestress",planestress, false);
181 :
182 6 : if (AMREX_SPACEDIM==2 && planestress)
183 5 : value.Define(mu,lambda*(1.0 - lambda/(2.*mu + lambda)));
184 : else
185 2 : value.Define(mu,lambda);
186 28 : }
187 :
188 : #define OP_CLASS Isotropic
189 : #define OP_VARS X(ddw)
190 : #include "Model/Solid/InClassOperators.H"
191 : };
192 : #include "Model/Solid/ExtClassOperators.H"
193 :
194 :
195 :
196 : }
197 : }
198 : }
199 :
200 : #endif
|