11#include <dune/common/float_cmp.hh>
25 "These are hard-coded function you should use the material library Ikarus::LinearElasticity")]]
inline Eigen::
30 C(0, 0) = C(1, 1) = 1;
31 C(0, 1) = C(1, 0) = nu;
32 C(2, 2) = (1 - nu) / 2;
33 C *= E / (1 - nu * nu);
46 "These are hard-coded function you should use the material library: Ikarus::LinearElasticity")]]
inline Eigen::
49 Eigen::Matrix<double, 6, 6> C;
51 C(0, 0) = C(1, 1) = C(2, 2) = 1 - nu;
52 C(0, 1) = C(1, 0) = C(2, 0) = C(0, 2) = C(1, 2) = C(2, 1) = nu;
53 C(3, 3) = C(4, 4) = C(5, 5) = (1 - 2 * nu) / 2;
54 C *= E / ((1 + nu) * (1 - 2 * nu));
64 template <
typename LocalView,
bool useRef = false>
68 static constexpr int worlddim = GridEntity::Geometry::coorddimension;
71 static constexpr int mydim = GridEntity::mydimension;
74 static constexpr int dimension = GridEntity::dimension;
77 template <
typename ScalarType =
double>
78 using VectorType = std::conditional_t<useRef, Eigen::Ref<Eigen::VectorX<ScalarType>>, Eigen::VectorX<ScalarType>&>;
81 template <
typename ScalarType =
double>
82 using MatrixType = std::conditional_t<useRef, Eigen::Ref<Eigen::MatrixX<ScalarType>>, Eigen::MatrixX<ScalarType>&>;
127 template <
typename MaterialParameter>
129 MaterialParameter, YoungsModulusAndBulkModulus> or std::is_same_v<MaterialParameter,
130 YoungsModulusAndLamesFirstParameter> or std::
131 is_same_v<MaterialParameter, BulkModulusAndLamesFirstParameter> or std::is_same_v<
132 MaterialParameter, LamesFirstParameterAndShearModulus> or std::is_same_v<MaterialParameter,
133 YoungsModulusAndShearModulus>;
140 template <
typename ValuePair>
145 YoungsModulusAndLamesFirstParameter> and !std::is_same_v<ValuePair, BulkModulusAndLamesFirstParameter> and !std::is_same_v<ValuePair, LamesFirstParameterAndShearModulus>) {
146 if constexpr (std::is_same_v<ValuePair, YoungsModulusAndPoissonsRatio>) {
147 const auto& E = vp.emodul;
148 const auto& nu = vp.nu;
149 return Dune::FloatCmp::eq(nu, 0.5) ? std::numeric_limits<double>::infinity()
150 : E * nu / ((1.0 + nu) * (1.0 - 2.0 * nu));
151 }
else if constexpr (std::is_same_v<ValuePair, YoungsModulusAndShearModulus>) {
152 const auto& E = vp.emodul;
153 const auto& mu = vp.mu;
154 return mu * (E - 2.0 * mu) / (3.0 * mu - E);
155 }
else if constexpr (std::is_same_v<ValuePair, YoungsModulusAndBulkModulus>) {
156 const auto& E = vp.emodul;
157 const auto& K = vp.K;
158 return 3.0 * K * (3.0 * K - E) / (9.0 * K - E);
160 assert(
false &&
"Your LameParameter request is not implemented");
166 if constexpr (std::is_same_v<ValuePair, YoungsModulusAndPoissonsRatio>) {
167 const auto& E = vp.emodul;
168 const auto& nu = vp.nu;
169 return E / (3.0 * (1.0 - 2.0 * nu));
170 }
else if constexpr (std::is_same_v<ValuePair, YoungsModulusAndShearModulus>) {
171 const auto& E = vp.emodul;
172 const auto& mu = vp.mu;
173 return E * mu / (3.0 * (3.0 * mu - E));
174 }
else if constexpr (std::is_same_v<ValuePair, YoungsModulusAndLamesFirstParameter>) {
175 const auto& E = vp.emodul;
176 const auto& lambda = vp.lambda;
177 return (E + 3.0 * lambda + calcR(vp)) / 6.0;
178 }
else if constexpr (std::is_same_v<ValuePair, LamesFirstParameterAndShearModulus>) {
179 const auto& lambda = vp.lambda;
180 const auto& mu = vp.mu;
181 return lambda + 2.0 * mu / 3.0;
183 assert(
false &&
"Your LameParameter request is not implemented");
190 if constexpr (std::is_same_v<ValuePair, YoungsModulusAndPoissonsRatio>) {
191 const auto& E = vp.emodul;
192 const auto& nu = vp.nu;
193 return E / (2.0 * (1.0 + nu));
194 }
else if constexpr (std::is_same_v<ValuePair, YoungsModulusAndBulkModulus>) {
195 const auto& E = vp.emodul;
196 const auto& K = vp.K;
197 return 3.0 * K * E / (9.0 * K - E);
198 }
else if constexpr (std::is_same_v<ValuePair, YoungsModulusAndLamesFirstParameter>) {
199 const auto& E = vp.emodul;
200 const auto& lambda = vp.lambda;
201 return (E - 3.0 * lambda + calcR(vp)) / 4.0;
202 }
else if constexpr (std::is_same_v<ValuePair, BulkModulusAndLamesFirstParameter>) {
203 const auto& K = vp.K;
204 const auto& lambda = vp.lambda;
205 return 3.0 * (K - lambda) / 2.0;
207 assert(
false &&
"Your LameParameter request is not implemented");
211 if constexpr (std::is_same_v<ValuePair, YoungsModulusAndPoissonsRatio>) {
212 const auto& E = vp.emodul;
213 const auto& nu = vp.nu;
214 return E * (1.0 - nu) / ((1.0 + nu) * (1.0 - 2.0 * nu));
215 }
else if constexpr (std::is_same_v<ValuePair, YoungsModulusAndShearModulus>) {
216 const auto& E = vp.emodul;
217 const auto& mu = vp.mu;
218 return mu * (4.0 * mu - E) / (3.0 * mu - E);
219 }
else if constexpr (std::is_same_v<ValuePair, YoungsModulusAndBulkModulus>) {
220 const auto& E = vp.emodul;
221 const auto& K = vp.K;
222 return 3.0 * K * (3.0 * K + E) / (9.0 * K - E);
223 }
else if constexpr (std::is_same_v<ValuePair, YoungsModulusAndLamesFirstParameter>) {
224 const auto& E = vp.emodul;
225 const auto& lambda = vp.lambda;
226 return (E - lambda + calcR(vp)) / 2.0;
227 }
else if constexpr (std::is_same_v<ValuePair, BulkModulusAndLamesFirstParameter>) {
228 const auto& K = vp.K;
229 const auto& lambda = vp.lambda;
230 return 3.0 * K - 2.0 * lambda;
231 }
else if constexpr (std::is_same_v<ValuePair, LamesFirstParameterAndShearModulus>) {
232 const auto& lambda = vp.lambda;
233 const auto& mu = vp.mu;
234 return lambda + 2.0 * mu;
236 assert(
false &&
"Your LameParameter request is not implemented");
239 constexpr double toPoissonsRatio()
requires(!std::is_same_v<ValuePair, YoungsModulusAndPoissonsRatio>) {
240 if constexpr (std::is_same_v<ValuePair, YoungsModulusAndShearModulus>) {
241 const auto& E = vp.emodul;
242 const auto& mu = vp.mu;
243 return E / (2.0 * mu) - 1.0;
244 }
else if constexpr (std::is_same_v<ValuePair, YoungsModulusAndBulkModulus>) {
245 const auto& E = vp.emodul;
246 const auto& K = vp.K;
247 return (3.0 * K - E) / (6.0 * K);
248 }
else if constexpr (std::is_same_v<ValuePair, YoungsModulusAndLamesFirstParameter>) {
249 const auto& E = vp.emodul;
250 const auto& lambda = vp.lambda;
251 return 2.0 * lambda / (E + lambda + calcR(vp));
252 }
else if constexpr (std::is_same_v<ValuePair, BulkModulusAndLamesFirstParameter>) {
253 const auto& K = vp.K;
254 const auto& lambda = vp.lambda;
255 return lambda / (3 * K - lambda);
256 }
else if constexpr (std::is_same_v<ValuePair, LamesFirstParameterAndShearModulus>) {
257 const auto& lambda = vp.lambda;
258 const auto& mu = vp.mu;
259 return lambda / (2.0 * (lambda + mu));
261 assert(
false &&
"Your LameParameter request is not implemented");
267 YoungsModulusAndPoissonsRatio> and !std::is_same_v<ValuePair, YoungsModulusAndShearModulus> and !std::is_same_v<ValuePair, YoungsModulusAndBulkModulus> and !std::is_same_v<ValuePair, YoungsModulusAndLamesFirstParameter>) {
268 if constexpr (std::is_same_v<ValuePair, BulkModulusAndLamesFirstParameter>) {
269 return 9.0 * vp.K * (vp.K - vp.lambda) / (3.0 * vp.K - vp.lambda);
270 }
else if constexpr (std::is_same_v<ValuePair, LamesFirstParameterAndShearModulus>) {
271 const auto& lambda = vp.lambda;
272 const auto& mu = vp.mu;
273 return mu * (3.0 * lambda + 2.0 * mu) / (lambda + mu);
275 assert(
false &&
"Your LameParameter request is not implemented");
292 ConvertLameConstants(
const ValuePair& p_vp) : vp(p_vp) {}
294 double calcR(
const YoungsModulusAndLamesFirstParameter& vp_) {
295 const auto& E = vp_.emodul;
296 const auto& lambda = vp_.lambda;
297 return std::sqrt(E * E + 9 * lambda * lambda + 2 * E * lambda);
Definition: simpleassemblers.hh:21
auto toLamesFirstParameterAndShearModulus(const YoungsModulusAndPoissonsRatio &matParameter)
Converts Young's modulus and Poisson's ratio to Lame's first parameter and shear modulus.
Definition: physicshelper.hh:328
auto toYoungsModulusAndPoissonsRatio(const LamesFirstParameterAndShearModulus &matParameter)
Converts Lame's first parameter and shear modulus to Young's modulus and Poisson's ratio.
Definition: physicshelper.hh:341
Eigen::Matrix3d planeStressLinearElasticMaterialTangent(double E, double nu)
Computes the plane stress linear elastic material tangent matrix.
Definition: physicshelper.hh:27
ConvertLameConstants< YoungsModulusAndPoissonsRatio > convertLameConstants(const YoungsModulusAndPoissonsRatio &p_vp)
Definition: physicshelper.hh:301
Eigen::Matrix< double, 6, 6 > linearElasticMaterialTangent3D(double E, double nu)
Computes the 3D linear elastic material tangent matrix.
Definition: physicshelper.hh:48
Traits for handling local views.see https://en.wikipedia.org/wiki/Lam%C3%A9_parameters.
Definition: physicshelper.hh:65
std::conditional_t< useRef, Eigen::Ref< Eigen::MatrixX< ScalarType > >, Eigen::MatrixX< ScalarType > & > MatrixType
Type of the stiffness matrix.
Definition: physicshelper.hh:82
static constexpr int worlddim
Dimension of the world space.
Definition: physicshelper.hh:68
std::conditional_t< useRef, Eigen::Ref< Eigen::VectorX< ScalarType > >, Eigen::VectorX< ScalarType > & > VectorType
Type of the internal forces.
Definition: physicshelper.hh:78
typename LocalView::Element GridEntity
Definition: physicshelper.hh:66
static constexpr int dimension
Dimension of the grid.
Definition: physicshelper.hh:74
static constexpr int mydim
Dimension of the geometry.
Definition: physicshelper.hh:71
Structure representing Young's modulus and shear modulus.
Definition: physicshelper.hh:87
double emodul
Definition: physicshelper.hh:88
double nu
Definition: physicshelper.hh:89
Structure representing Young's modulus and bulk modulus.
Definition: physicshelper.hh:93
double mu
Definition: physicshelper.hh:95
double emodul
Definition: physicshelper.hh:94
Structure representing Young's modulus and Lame's first parameter.
Definition: physicshelper.hh:99
double emodul
Definition: physicshelper.hh:100
double K
Definition: physicshelper.hh:101
Structure representing bulk modulus and Lame's first parameter.
Definition: physicshelper.hh:105
double lambda
Definition: physicshelper.hh:107
double emodul
Definition: physicshelper.hh:106
Structure representing Lame's first parameter and shear modulus.
Definition: physicshelper.hh:111
double K
Definition: physicshelper.hh:112
double lambda
Definition: physicshelper.hh:113
Definition: physicshelper.hh:117
double lambda
Definition: physicshelper.hh:118
double mu
Definition: physicshelper.hh:119
Conversion utility for Lame's constants.
Definition: physicshelper.hh:141
constexpr double toBulkModulus()
Definition: physicshelper.hh:163
constexpr double toYoungsModulus()
Definition: physicshelper.hh:264
constexpr double toPoissonsRatio()
Definition: physicshelper.hh:239
friend ConvertLameConstants< YoungsModulusAndPoissonsRatio > convertLameConstants(const YoungsModulusAndPoissonsRatio &p_vp)
Definition: physicshelper.hh:301
constexpr double toLamesFirstParameter()
Definition: physicshelper.hh:142
constexpr double toPWaveModulus()
Definition: physicshelper.hh:210
constexpr double toShearModulus()
Definition: physicshelper.hh:186
Concept for checking if a type is a valid material parameter tuple.
Definition: physicshelper.hh:128