OGS
LinearElasticOrthotropic.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <limits>
7
8#include "MechanicsBase.h"
10
11namespace MaterialLib
12{
13namespace Solids
14{
15template <int DisplacementDim>
16class LinearElasticOrthotropic : public MechanicsBase<DisplacementDim>
17{
18public:
20 {
22 constexpr double E(int const i) const
23 {
24 assert(i == 1 || i == 2 || i == 3);
25 return youngs_moduli[i - 1];
26 }
27
29 constexpr double G(int const i, int const j) const
30 {
31 assert(i == 1 || i == 2 || i == 3);
32 assert(j == 1 || j == 2 || j == 3);
33 if (i == 1 && j == 2)
34 {
35 return shear_moduli[0];
36 }
37 if (i == 2 && j == 3)
38 {
39 return shear_moduli[1];
40 }
41 if (i == 1 && j == 3)
42 {
43 return shear_moduli[2];
44 }
45
46 return std::numeric_limits<double>::quiet_NaN();
47 }
48
50 constexpr double nu(int const i, int const j) const
51 {
52 assert(i == 1 || i == 2 || i == 3);
53 assert(j == 1 || j == 2 || j == 3);
54 if (i == 1 && j == 2)
55 {
56 return poissons_ratios[0];
57 }
58 if (i == 2 && j == 3)
59 {
60 return poissons_ratios[1];
61 }
62 if (i == 1 && j == 3)
63 {
64 return poissons_ratios[2];
65 }
66
67 // The next set is scaled by Youngs modulus s.t.
68 // nu_ji = nu_ij * E_j/E_i.
69 if (i == 2 && j == 1)
70 {
71 return poissons_ratios[0] * E(i) / E(j);
72 }
73 if (i == 3 && j == 2)
74 {
75 return poissons_ratios[1] * E(i) / E(j);
76 }
77 if (i == 3 && j == 1)
78 {
79 return poissons_ratios[2] * E(i) / E(j);
80 }
81
82 return std::numeric_limits<double>::quiet_NaN();
83 }
84
85 std::array<double, 3> youngs_moduli;
86 std::array<double, 3> shear_moduli;
87 std::array<double, 3> poissons_ratios;
88 };
89
92 {
95
96 MaterialProperties(P const& youngs_moduli_, P const& shear_moduli_,
97 P const& poissons_ratios_)
98 : youngs_moduli(youngs_moduli_),
99 shear_moduli(shear_moduli_),
100 poissons_ratios(poissons_ratios_)
101 {
102 }
103
105 double const t, ParameterLib::SpatialPosition const& x) const
106 {
107 auto const E = youngs_moduli(t, x);
108 auto const G = shear_moduli(t, x);
109 auto const nu = poissons_ratios(t, x);
110 return {
111 {E[0], E[1], E[2]}, {G[0], G[1], G[2]}, {nu[0], nu[1], nu[2]}};
112 }
113
116 P const& poissons_ratios; // Stored as nu_12, nu_23, nu_13
117 };
118
119public:
120 static int const KelvinVectorSize =
126
128 MaterialProperties material_properties,
129 std::optional<ParameterLib::CoordinateSystem> const&
130 local_coordinate_system)
131 : _mp(std::move(material_properties)),
132 _local_coordinate_system(local_coordinate_system)
133 {
134 }
135
137 double const /*t*/,
139 double const /*dt*/,
140 KelvinVector const& eps,
141 KelvinVector const& sigma,
143 MaterialStateVariables const& /* material_state_variables */)
144 const override
145 {
146 return eps.dot(sigma) / 2;
147 }
148
149 std::optional<
150 std::tuple<typename MechanicsBase<DisplacementDim>::KelvinVector,
151 std::unique_ptr<typename MechanicsBase<
152 DisplacementDim>::MaterialStateVariables>,
155 MaterialPropertyLib::VariableArray const& variable_array_prev,
156 MaterialPropertyLib::VariableArray const& variable_array,
157 double const t, ParameterLib::SpatialPosition const& x,
158 double const /*dt*/,
160 material_state_variables) const override;
161
162 KelvinMatrix getElasticTensor(double const t,
164 double const T) const;
165
167
168 double getBulkModulus(double const t,
170 KelvinMatrix const* const /*C*/) const override
171 {
172 auto const& mp = _mp.evaluate(t, x);
173 auto const E = [&mp](int const i) { return mp.E(i); };
174 auto const nu = [&mp](int const i, int const j) { return mp.nu(i, j); };
175 // corresponds to 1/(I:S:I) --> Reuss bound.
176 // Voigt bound would proceed as (I:C:I)/9. Use MFront model if you
177 // prefer that.
178 return E(1) * E(2) * E(3) /
179 (E(1) * E(2) + E(1) * E(3) * (1 - 2 * nu(2, 3)) +
180 E(2) * E(3) * (1 - 2 * nu(1, 2) * nu(1, 3)));
181 }
182
183protected:
185 std::optional<ParameterLib::CoordinateSystem> const&
187};
188
189extern template class LinearElasticOrthotropic<2>;
190extern template class LinearElasticOrthotropic<3>;
191
192} // namespace Solids
193} // namespace MaterialLib
double getBulkModulus(double const t, ParameterLib::SpatialPosition const &x, KelvinMatrix const *const) const override
LinearElasticOrthotropic(MaterialProperties material_properties, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system)
double computeFreeEnergyDensity(double const, ParameterLib::SpatialPosition const &, double const, KelvinVector const &eps, KelvinVector const &sigma, typename MechanicsBase< DisplacementDim >::MaterialStateVariables const &) const override
std::optional< ParameterLib::CoordinateSystem > const & _local_coordinate_system
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > KelvinMatrix
std::optional< std::tuple< typename MechanicsBase< DisplacementDim >::KelvinVector, std::unique_ptr< typename MechanicsBase< DisplacementDim >::MaterialStateVariables >, typename MechanicsBase< DisplacementDim >::KelvinMatrix > > integrateStress(MaterialPropertyLib::VariableArray const &variable_array_prev, MaterialPropertyLib::VariableArray const &variable_array, double const t, ParameterLib::SpatialPosition const &x, double const, typename MechanicsBase< DisplacementDim >::MaterialStateVariables const &material_state_variables) const override
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > KelvinVector
KelvinMatrix getElasticTensor(double const t, ParameterLib::SpatialPosition const &x, double const T) const
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType
constexpr double G(int const i, int const j) const
Shear moduli for 1-based indices.
constexpr double E(int const i) const
Youngs moduli for 1-based indices.
constexpr double nu(int const i, int const j) const
Poisson's ratios for 1-based indices.
EvaluatedMaterialProperties evaluate(double const t, ParameterLib::SpatialPosition const &x) const
MaterialProperties(P const &youngs_moduli_, P const &shear_moduli_, P const &poissons_ratios_)
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > KelvinMatrix