OGS
SmallDeformationFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <limits>
14#include <memory>
15#include <vector>
16
21#include "NumLib/DOF/LocalDOF.h"
32
33namespace ProcessLib
34{
35namespace SmallDeformation
36{
37namespace MPL = MaterialPropertyLib;
38
39template <typename BMatricesType, typename ShapeMatricesType,
40 int DisplacementDim>
42{
44 typename ShapeMatricesType::NodalRowVectorType N;
45 typename ShapeMatricesType::GlobalDimNodalMatrixType dNdx;
46
48};
49
52template <typename ShapeMatrixType>
54{
55 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N;
56};
57
58template <typename ShapeFunction, int DisplacementDim>
60 : public SmallDeformationLocalAssemblerInterface<DisplacementDim>
61{
62public:
69
75
79 using IpData =
81
82 static constexpr auto& N_u_op = MathLib::eigenBlockMatrixView<
83 DisplacementDim, typename ShapeMatricesType::NodalRowVectorType>;
84
86 delete;
88
90 MeshLib::Element const& e,
91 std::size_t const /*local_matrix_size*/,
92 NumLib::GenericIntegrationMethod const& integration_method,
93 bool const is_axially_symmetric,
96 e, integration_method, is_axially_symmetric, process_data)
97 {
98 unsigned const n_integration_points =
100
101 _ip_data.resize(n_integration_points);
102 _secondary_data.N.resize(n_integration_points);
103
104 auto const shape_matrices =
106 DisplacementDim>(
107 e, is_axially_symmetric, this->integration_method_);
108
109 for (unsigned ip = 0; ip < n_integration_points; ip++)
110 {
111 auto& ip_data = _ip_data[ip];
112 auto const& sm = shape_matrices[ip];
113 _ip_data[ip].integration_weight =
115 sm.integralMeasure * sm.detJ;
116
117 ip_data.N = sm.N;
118 ip_data.dNdx = sm.dNdx;
119
120 _secondary_data.N[ip] = shape_matrices[ip].N;
121 }
122 }
123
124 void initializeConcrete() override
125 {
126 unsigned const n_integration_points =
128 for (unsigned ip = 0; ip < n_integration_points; ip++)
129 {
130 auto const& ip_data = _ip_data[ip];
131
132 ParameterLib::SpatialPosition const x_position{
133 std::nullopt, this->element_.getID(), ip,
135 NumLib::interpolateCoordinates<ShapeFunction,
137 this->element_, ip_data.N))};
138
140 if (this->process_data_.initial_stress != nullptr)
141 {
142 this->current_states_[ip].stress_data.sigma.noalias() =
144 DisplacementDim>((*this->process_data_.initial_stress)(
145 std::numeric_limits<
146 double>::quiet_NaN() /* time independent */,
147 x_position));
148 }
149
150 double const t = 0; // TODO (naumov) pass t from top
151 auto& material_state = this->material_states_[ip];
152 this->solid_material_.initializeInternalStateVariables(
153 t, x_position, *material_state.material_state_variables);
154
155 material_state.pushBackState();
156 }
157
158 for (unsigned ip = 0; ip < n_integration_points; ip++)
159 {
160 this->prev_states_[ip] = this->current_states_[ip];
161 }
162 }
163
166 Eigen::Ref<Eigen::VectorXd const> const& u,
167 Eigen::Ref<Eigen::VectorXd const> const& u_prev,
168 ParameterLib::SpatialPosition const& x_position, double const t,
169 double const dt,
171 ip_data,
173 CS,
174 MaterialPropertyLib::Medium const& medium,
176 current_state,
178 prev_state,
181 output_data) const
182 {
183 auto const& N = ip_data.N;
184 auto const& dNdx = ip_data.dNdx;
185 auto const x_coord =
186 NumLib::interpolateXCoordinate<ShapeFunction, ShapeMatricesType>(
187 this->element_, N);
188 auto const B =
189 LinearBMatrix::computeBMatrix<DisplacementDim,
190 ShapeFunction::NPOINTS,
192 dNdx, N, x_coord, this->is_axially_symmetric_);
193
194 double const T_ref =
195 this->process_data_.reference_temperature
196 ? (*this->process_data_.reference_temperature)(t, x_position)[0]
197 : std::numeric_limits<double>::quiet_NaN();
198
200 models{this->process_data_, this->solid_material_};
202 tmp;
204
205 CS.eval(models, t, dt, x_position, //
206 medium, //
207 T_ref, B * u, B * u_prev, //
208 current_state, prev_state, material_state, tmp, output_data,
209 CD);
210
211 return CD;
212 }
213
214 void assemble(double const /*t*/, double const /*dt*/,
215 std::vector<double> const& /*local_x*/,
216 std::vector<double> const& /*local_x_prev*/,
217 std::vector<double>& /*local_M_data*/,
218 std::vector<double>& /*local_K_data*/,
219 std::vector<double>& /*local_b_data*/) override
220 {
221 OGS_FATAL(
222 "SmallDeformationLocalAssembler: assembly without jacobian is not "
223 "implemented.");
224 }
225
226 void assembleWithJacobian(double const t, double const dt,
227 std::vector<double> const& local_x,
228 std::vector<double> const& local_x_prev,
229 std::vector<double>& /*local_M_data*/,
230 std::vector<double>& /*local_K_data*/,
231 std::vector<double>& local_b_data,
232 std::vector<double>& local_Jac_data) override
233 {
234 auto const local_matrix_size = local_x.size();
235
236 auto local_Jac = MathLib::createZeroedMatrix<StiffnessMatrixType>(
237 local_Jac_data, local_matrix_size, local_matrix_size);
238
239 auto local_b = MathLib::createZeroedVector<NodalDisplacementVectorType>(
240 local_b_data, local_matrix_size);
241
242 auto [u] = localDOF(local_x);
243 auto [u_prev] = localDOF(local_x_prev);
244
245 unsigned const n_integration_points =
247
249 x_position.setElementID(this->element_.getID());
250
253 auto const& medium =
254 *this->process_data_.media_map.getMedium(this->element_.getID());
255
256 for (unsigned ip = 0; ip < n_integration_points; ip++)
257 {
258 x_position.setIntegrationPoint(ip);
259 auto const& w = _ip_data[ip].integration_weight;
260 auto const& N = _ip_data[ip].N;
261 auto const& dNdx = _ip_data[ip].dNdx;
262
263 auto const x_coord =
264 NumLib::interpolateXCoordinate<ShapeFunction,
266 this->element_, N);
267 auto const B = LinearBMatrix::computeBMatrix<
268 DisplacementDim, ShapeFunction::NPOINTS,
270 dNdx, N, x_coord, this->is_axially_symmetric_);
271
272 auto const CD = updateConstitutiveRelations(
273 u, u_prev, x_position, t, dt, _ip_data[ip],
274 constitutive_setting, medium, this->current_states_[ip],
275 this->prev_states_[ip], this->material_states_[ip],
276 this->output_data_[ip]);
277
278 auto const& sigma = this->current_states_[ip].stress_data.sigma;
279 auto const& b = *CD.volumetric_body_force;
280 auto const& C = CD.s_mech_data.stiffness_tensor;
281
282 local_b.noalias() -=
283 (B.transpose() * sigma - N_u_op(N).transpose() * b) * w;
284 local_Jac.noalias() += B.transpose() * C * B * w;
285 }
286 }
287
288 void postTimestepConcrete(Eigen::VectorXd const& local_x,
289 Eigen::VectorXd const& local_x_prev,
290 double const t, double const dt,
291 int const /*process_id*/) override
292 {
293 unsigned const n_integration_points =
295
297 x_position.setElementID(this->element_.getID());
298
301
302 auto const& medium =
303 *this->process_data_.media_map.getMedium(this->element_.getID());
304
305 for (unsigned ip = 0; ip < n_integration_points; ip++)
306 {
307 x_position.setIntegrationPoint(ip);
308
310 local_x, local_x_prev, x_position, t, dt, _ip_data[ip],
311 constitutive_setting, medium, this->current_states_[ip],
312 this->prev_states_[ip], this->material_states_[ip],
313 this->output_data_[ip]);
314
315 this->material_states_[ip].pushBackState();
316 }
317
318 for (unsigned ip = 0; ip < n_integration_points; ip++)
319 {
320 this->prev_states_[ip] = this->current_states_[ip];
321 }
322 }
323
324 std::vector<double> const& getMaterialForces(
325 std::vector<double> const& local_x,
326 std::vector<double>& nodal_values) override
327 {
329 DisplacementDim, ShapeFunction, ShapeMatricesType,
332 GradientMatrixType>(local_x, nodal_values,
334 this->current_states_, this->output_data_,
335 this->element_, this->is_axially_symmetric_);
336 }
337
338 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
339 const unsigned integration_point) const override
340 {
341 auto const& N = _secondary_data.N[integration_point];
342
343 // assumes N is stored contiguously in memory
344 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
345 }
346
347private:
348 static constexpr auto localDOF(std::vector<double> const& x)
349 {
350 return NumLib::localDOF<
352 }
353
354private:
355 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
356
358};
359
360} // namespace SmallDeformation
361} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
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
MatrixType< _number_of_dof, _number_of_dof > StiffnessMatrixType
VectorType< _number_of_dof > NodalForceVectorType
Rhs residual.
VectorTypeFixedSize< MathLib::VectorizedTensor::size(DisplacementDim)> GradientVectorType
MatrixType< MathLib::VectorizedTensor::size(DisplacementDim), _number_of_dof > GradientMatrixType
typename BMatricesType::NodalForceVectorType NodalForceVectorType
typename ShapeMatricesType::NodalVectorType NodalVectorType
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
ConstitutiveRelations::ConstitutiveData< DisplacementDim > updateConstitutiveRelations(Eigen::Ref< Eigen::VectorXd const > const &u, Eigen::Ref< Eigen::VectorXd const > const &u_prev, ParameterLib::SpatialPosition const &x_position, double const t, double const dt, IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim > &ip_data, typename ConstitutiveRelations::ConstitutiveSetting< DisplacementDim > &CS, MaterialPropertyLib::Medium const &medium, typename ConstitutiveRelations::StatefulData< DisplacementDim > &current_state, typename ConstitutiveRelations::StatefulDataPrev< DisplacementDim > const &prev_state, MaterialStateData< DisplacementDim > &material_state, typename ConstitutiveRelations::OutputData< DisplacementDim > &output_data) const
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
typename BMatricesType::NodalForceVectorType NodalDisplacementVectorType
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
void postTimestepConcrete(Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, double const t, double const dt, int const) override
typename ShapeMatricesType::NodalMatrixType NodalMatrixType
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &, std::vector< double > &, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
typename GMatricesType::GradientMatrixType GradientMatrixType
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
std::vector< double > const & getMaterialForces(std::vector< double > const &local_x, std::vector< double > &nodal_values) override
typename BMatricesType::StiffnessMatrixType StiffnessMatrixType
SmallDeformationLocalAssembler(SmallDeformationLocalAssembler &&)=delete
static constexpr auto localDOF(std::vector< double > const &x)
typename GMatricesType::GradientVectorType GradientVectorType
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
SmallDeformationLocalAssembler(MeshLib::Element const &e, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, SmallDeformationProcessData< DisplacementDim > &process_data)
SmallDeformationLocalAssembler(SmallDeformationLocalAssembler const &)=delete
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)
constexpr Eigen::CwiseNullaryOp< EigenBlockMatrixViewFunctor< D, M >, typename EigenBlockMatrixViewFunctor< D, M >::Matrix > eigenBlockMatrixView(const Eigen::MatrixBase< M > &matrix)
double interpolateXCoordinate(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
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)
auto localDOF(ElementDOFVector const &x)
Definition LocalDOF.h:64
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
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< double > const & getMaterialForces(std::vector< double > const &local_x, std::vector< double > &nodal_values, IntegrationMethod const &_integration_method, IPData const &_ip_data, StressData const &stress_data, OutputData const &output_data, MeshLib::Element const &element, bool const is_axially_symmetric)
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N
Data that is needed for the equation system assembly.
void eval(ConstitutiveModels< DisplacementDim > &models, double const t, double const dt, ParameterLib::SpatialPosition const &x_position, MaterialPropertyLib::Medium const &medium, double const T_ref, KelvinVector< DisplacementDim > const &eps, KelvinVector< DisplacementDim > const &eps_prev, StatefulData< DisplacementDim > &state, StatefulDataPrev< DisplacementDim > const &prev_state, MaterialStateData< DisplacementDim > &mat_state, ConstitutiveTempData< DisplacementDim > &tmp, OutputData< DisplacementDim > &out, ConstitutiveData< DisplacementDim > &cd) const
Evaluate the constitutive setting.
Data that is needed for output purposes, but not directly for the assembly.
Data whose state must be tracked by the process.
ShapeMatricesType::GlobalDimNodalMatrixType dNdx
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N
ConstitutiveRelations::ConstitutiveSetting< DisplacementDim > constitutive_setting
std::vector< typename ConstitutiveRelations::StatefulDataPrev< DisplacementDim > > prev_states_
std::vector< typename ConstitutiveRelations::StatefulData< DisplacementDim > > current_states_
std::vector< typename ConstitutiveRelations::OutputData< DisplacementDim > > output_data_
MaterialLib::Solids::MechanicsBase< DisplacementDim > const & solid_material_
std::vector< MaterialStateData< DisplacementDim > > material_states_