OGS
SmallDeformationLocalAssemblerMatrix-impl.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <Eigen/Core>
14#include <valarray>
15#include <vector>
16
28#include "SecondaryData.h"
31
32namespace ProcessLib
33{
34namespace LIE
35{
36namespace SmallDeformation
37{
38template <typename ShapeFunction, int DisplacementDim>
41 MeshLib::Element const& e,
42 std::size_t const /*local_matrix_size*/,
43 NumLib::GenericIntegrationMethod const& integration_method,
44 bool const is_axially_symmetric,
46 : _process_data(process_data),
47 _integration_method(integration_method),
48 _element(e),
49 _is_axially_symmetric(is_axially_symmetric)
50{
51 unsigned const n_integration_points =
53
54 _ip_data.reserve(n_integration_points);
55 _secondary_data.N.resize(n_integration_points);
56
57 auto const shape_matrices =
59 DisplacementDim>(e, is_axially_symmetric,
61
63 _process_data.solid_materials, _process_data.material_ids, e.getID());
64
65 for (unsigned ip = 0; ip < n_integration_points; ip++)
66 {
67 _ip_data.emplace_back(solid_material);
68 auto& ip_data = _ip_data[ip];
69 auto const& sm = shape_matrices[ip];
70 ip_data.N = sm.N;
71 ip_data.dNdx = sm.dNdx;
72 ip_data.integration_weight =
74 sm.integralMeasure * sm.detJ;
75
76 // Initialize current time step values
77 static const int kelvin_vector_size =
79 ip_data._sigma.setZero(kelvin_vector_size);
80 ip_data._eps.setZero(kelvin_vector_size);
81
82 // Previous time step values are not initialized and are set later.
83 ip_data._sigma_prev.resize(kelvin_vector_size);
84 ip_data._eps_prev.resize(kelvin_vector_size);
85
86 ip_data._C.resize(kelvin_vector_size, kelvin_vector_size);
87
88 _secondary_data.N[ip] = sm.N;
89 }
90}
91
92template <typename ShapeFunction, int DisplacementDim>
94 assembleWithJacobian(double const t, double const dt,
95 std::vector<double> const& local_x,
96 std::vector<double> const& /*local_x_prev*/,
97 std::vector<double>& /*local_M_data*/,
98 std::vector<double>& /*local_K_data*/,
99 std::vector<double>& local_b_data,
100 std::vector<double>& local_Jac_data)
101{
102 assert(_element.getDimension() == DisplacementDim);
103
104 auto const local_matrix_size = local_x.size();
105
106 auto local_Jac = MathLib::createZeroedMatrix<StiffnessMatrixType>(
107 local_Jac_data, local_matrix_size, local_matrix_size);
108
109 auto local_b = MathLib::createZeroedVector<NodalDisplacementVectorType>(
110 local_b_data, local_matrix_size);
111
112 unsigned const n_integration_points =
113 _integration_method.getNumberOfPoints();
114
115 MPL::VariableArray variables;
116 MPL::VariableArray variables_prev;
118 x_position.setElementID(_element.getID());
119
120 for (unsigned ip = 0; ip < n_integration_points; ip++)
121 {
122 x_position.setIntegrationPoint(ip);
123 auto const& w = _ip_data[ip].integration_weight;
124
125 auto const& N = _ip_data[ip].N;
126 auto const& dNdx = _ip_data[ip].dNdx;
127 auto const x_coord =
128 NumLib::interpolateXCoordinate<ShapeFunction, ShapeMatricesType>(
129 _element, N);
130 auto const B =
131 LinearBMatrix::computeBMatrix<DisplacementDim,
132 ShapeFunction::NPOINTS,
134 dNdx, N, x_coord, _is_axially_symmetric);
135
136 auto const& eps_prev = _ip_data[ip]._eps_prev;
137 auto const& sigma_prev = _ip_data[ip]._sigma_prev;
138
139 auto& eps = _ip_data[ip]._eps;
140 auto& sigma = _ip_data[ip]._sigma;
141 auto& state = _ip_data[ip]._material_state_variables;
142
143 eps.noalias() =
144 B * Eigen::Map<typename BMatricesType::NodalForceVectorType const>(
145 local_x.data(), ShapeFunction::NPOINTS * DisplacementDim);
146
147 variables.mechanical_strain
149 eps);
150 variables_prev.stress
152 sigma_prev);
153 variables_prev.mechanical_strain
155 eps_prev);
156 variables_prev.temperature = _process_data._reference_temperature;
157
158 auto&& solution = _ip_data[ip]._solid_material.integrateStress(
159 variables_prev, variables, t, x_position, dt, *state);
160
161 if (!solution)
162 {
163 OGS_FATAL("Computation of local constitutive relation failed.");
164 }
165
167 std::tie(sigma, state, C) = std::move(*solution);
168
169 local_b.noalias() -= B.transpose() * sigma * w;
170 local_Jac.noalias() += B.transpose() * C * B * w;
171 }
172}
173
174template <typename ShapeFunction, int DisplacementDim>
177 double const /*t*/, Eigen::VectorXd const& /*local_x*/)
178{
179 // Compute average value per element
181 KV sigma_avg = KV::Zero();
182 auto const e_id = _element.getID();
183
184 unsigned const n_integration_points =
185 _integration_method.getNumberOfPoints();
186 for (unsigned ip = 0; ip < n_integration_points; ip++)
187 {
188 sigma_avg += _ip_data[ip]._sigma;
189 }
190 sigma_avg /= n_integration_points;
191
192 Eigen::Map<KV>(
193 &(*_process_data.element_stresses)[e_id * KV::RowsAtCompileTime]) =
195}
196
197template <typename ShapeFunction, int DisplacementDim>
198std::vector<double> const&
201 const double /*t*/,
202 std::vector<GlobalVector*> const& /*x*/,
203 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
204 std::vector<double>& cache) const
205{
206 return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
207 _ip_data, &IntegrationPointDataType::_sigma, cache);
208}
209template <typename ShapeFunction, int DisplacementDim>
210std::vector<double> const&
213 const double /*t*/,
214 std::vector<GlobalVector*> const& /*x*/,
215 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
216 std::vector<double>& cache) const
217{
218 return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
219 _ip_data, &IntegrationPointDataType::_eps, cache);
220}
221
222} // namespace SmallDeformation
223} // namespace LIE
224} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
std::variant< std::monostate, Eigen::Matrix< double, 4, 1 >, Eigen::Matrix< double, 6, 1 > > mechanical_strain
std::variant< std::monostate, Eigen::Matrix< double, 4, 1 >, Eigen::Matrix< double, 6, 1 > > stress
double getWeight() const
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
void computeSecondaryVariableConcreteWithVector(double const t, Eigen::VectorXd const &local_x) override
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
std::vector< double > const & getIntPtSigma(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
std::vector< double > const & getIntPtEpsilon(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
SmallDeformationLocalAssemblerMatrix(SmallDeformationLocalAssemblerMatrix const &)=delete
auto & selectSolidConstitutiveRelation(SolidMaterialsMap const &constitutive_relations, MeshLib::PropertyVector< int > const *const material_ids, std::size_t const element_id)
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
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
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
BMatrixType computeBMatrix(DNDX_Type const &dNdx, N_Type const &N, const double radius, const bool is_axially_symmetric)
Fills a B-matrix based on given shape function dN/dx values.
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N