OGS
LargeDeformationFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <limits>
14#include <memory>
15#include <vector>
16
23#include "NumLib/DOF/LocalDOF.h"
36
37namespace ProcessLib
38{
39namespace LargeDeformation
40{
41namespace MPL = MaterialPropertyLib;
42
45template <typename ShapeMatrixType>
47{
48 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N;
49};
50
51template <int DisplacementDim, typename ShapeMatricesType>
52Eigen::Matrix<double, MPL::tensorSize(DisplacementDim),
53 MPL::tensorSize(DisplacementDim)>
54computeSigmaGeom(Eigen::Matrix3d const& sigma_tensor)
55{
56 static constexpr auto& sigma_geom_op = MathLib::eigenBlockMatrixView<
57 DisplacementDim,
58 Eigen::Matrix<double, DisplacementDim, DisplacementDim>>;
59
60 using SigmaGeom = Eigen::Matrix<double, MPL::tensorSize(DisplacementDim),
61 MPL::tensorSize(DisplacementDim)>;
62 if constexpr (DisplacementDim == 2)
63 {
64 SigmaGeom sigma_geom = SigmaGeom::Zero(5, 5);
65 sigma_geom.template block<4, 4>(0, 0) =
66 sigma_geom_op(sigma_tensor.template block<2, 2>(0, 0).eval());
67 sigma_geom(4, 4) = sigma_tensor(2, 2);
68
69 return sigma_geom;
70 }
71
72 if constexpr (DisplacementDim == 3)
73 {
74 return sigma_geom_op(sigma_tensor);
75 }
76}
77
78template <typename ShapeFunction, int DisplacementDim>
80 : public LargeDeformationLocalAssemblerInterface<DisplacementDim>
81{
82public:
89
95
99 using IpData =
101
102 static constexpr auto& N_u_op = MathLib::eigenBlockMatrixView<
103 DisplacementDim, typename ShapeMatricesType::NodalRowVectorType>;
104
106 delete;
108
110 MeshLib::Element const& e,
111 std::size_t const /*local_matrix_size*/,
112 NumLib::GenericIntegrationMethod const& integration_method,
113 bool const is_axially_symmetric,
116 e, integration_method, is_axially_symmetric, process_data)
117 {
118 unsigned const n_integration_points =
120
121 _ip_data.resize(n_integration_points);
122 _secondary_data.N.resize(n_integration_points);
123
124 auto const shape_matrices =
126 DisplacementDim>(
127 e, is_axially_symmetric, this->integration_method_);
128
129 for (unsigned ip = 0; ip < n_integration_points; ip++)
130 {
131 auto& ip_data = _ip_data[ip];
132 auto const& sm = shape_matrices[ip];
133 _ip_data[ip].integration_weight =
135 sm.integralMeasure * sm.detJ;
136
137 ip_data.N = sm.N;
138 ip_data.dNdx = sm.dNdx;
139
140 _secondary_data.N[ip] = shape_matrices[ip].N;
141 }
142 }
143
144 void initializeConcrete() override
145 {
146 unsigned const n_integration_points =
148 for (unsigned ip = 0; ip < n_integration_points; ip++)
149 {
150 auto const& ip_data = _ip_data[ip];
151
152 ParameterLib::SpatialPosition const x_position{
153 std::nullopt, this->element_.getID(), ip,
155 NumLib::interpolateCoordinates<ShapeFunction,
157 this->element_, ip_data.N))};
158
160 if (this->process_data_.initial_stress != nullptr)
161 {
162 this->current_states_[ip].stress_data.sigma.noalias() =
164 DisplacementDim>((*this->process_data_.initial_stress)(
165 std::numeric_limits<
166 double>::quiet_NaN() /* time independent */,
167 x_position));
168 }
169
170 double const t = 0; // TODO (naumov) pass t from top
171 auto& material_state = this->material_states_[ip];
172 this->solid_material_.initializeInternalStateVariables(
173 t, x_position, *material_state.material_state_variables);
174
175 material_state.pushBackState();
176 this->prev_states_[ip] = this->current_states_[ip];
177 }
178 }
179
182 BMatrixType const&
183 B /*B_{linear}+B_{nonlinear}/2 for Green-Lagrange strain*/,
184 GradientMatrixType const& G, GradientVectorType const& grad_u,
185 Eigen::Ref<Eigen::VectorXd const> const& u,
186 Eigen::Ref<Eigen::VectorXd const> const& u_prev,
187 ParameterLib::SpatialPosition const& x_position, double const t,
188 double const dt,
190 CS,
191 MaterialPropertyLib::Medium const& medium,
193 current_state,
195 prev_state,
198 output_data) const
199 {
200 double const T_ref =
201 this->process_data_.reference_temperature
202 ? (*this->process_data_.reference_temperature)(t, x_position)[0]
203 : std::numeric_limits<double>::quiet_NaN();
204
206 models(this->process_data_, this->solid_material_);
208 tmp;
210
211 // Note: Here B=B_{linear}+B_{nonlinear}/2 For Green-Lagrange strain.
212 output_data.eps_data.eps = B * u;
213 output_data.deformation_gradient_data.deformation_gradient =
215 output_data.deformation_gradient_data.volume_ratio =
217 output_data.deformation_gradient_data.deformation_gradient);
218
219 CS.eval(
220 models, t, dt, x_position, //
221 medium, //
222 T_ref, //
223 output_data.deformation_gradient_data, //
225 current_state, prev_state, material_state, tmp, CD);
226
227 return CD;
228 }
229
230 void assemble(double const /*t*/, double const /*dt*/,
231 std::vector<double> const& /*local_x*/,
232 std::vector<double> const& /*local_x_prev*/,
233 std::vector<double>& /*local_M_data*/,
234 std::vector<double>& /*local_K_data*/,
235 std::vector<double>& /*local_b_data*/) override
236 {
237 OGS_FATAL(
238 "LargeDeformationLocalAssembler: assembly without jacobian is not "
239 "implemented.");
240 }
241
242 void assembleWithJacobian(double const t, double const dt,
243 std::vector<double> const& local_x,
244 std::vector<double> const& local_x_prev,
245 std::vector<double>& local_b_data,
246 std::vector<double>& local_Jac_data) override
247 {
248 auto const local_matrix_size = local_x.size();
249
251 local_Jac_data, local_matrix_size, local_matrix_size);
252
254 local_b_data, local_matrix_size);
255
256 auto [u] = localDOF(local_x);
257 auto [u_prev] = localDOF(local_x_prev);
258
259 unsigned const n_integration_points =
261
263 x_position.setElementID(this->element_.getID());
264
266 constitutive_setting;
267 auto const& medium =
268 *this->process_data_.media_map.getMedium(this->element_.getID());
269
270 for (unsigned ip = 0; ip < n_integration_points; ip++)
271 {
272 x_position.setIntegrationPoint(ip);
273 auto const& w = _ip_data[ip].integration_weight;
274 auto const& N = _ip_data[ip].N;
275 auto const& dNdx = _ip_data[ip].dNdx;
276
277 auto const x_coord =
278 NumLib::interpolateXCoordinate<ShapeFunction,
280 this->element_, N);
281
282 // For the 2D case the 33-component is needed (and the four entries
283 // of the non-symmetric matrix); In 3d there are nine entries.
284 GradientMatrixType G(DisplacementDim * DisplacementDim +
285 (DisplacementDim == 2 ? 1 : 0),
286 ShapeFunction::NPOINTS * DisplacementDim);
287 Deformation::computeGMatrix<DisplacementDim,
288 ShapeFunction::NPOINTS>(
289 dNdx, G, this->is_axially_symmetric_, N, x_coord);
290
291 GradientVectorType const grad_u = G * u;
292
293 auto const B_0 = LinearBMatrix::computeBMatrix<
294 DisplacementDim, ShapeFunction::NPOINTS,
296 dNdx, N, x_coord, this->is_axially_symmetric_);
297
298 auto const B_N = NonLinearBMatrix::computeBMatrix<
299 DisplacementDim, ShapeFunction::NPOINTS,
301 dNdx, N, grad_u, x_coord, this->is_axially_symmetric_);
302
303 auto const CD = updateConstitutiveRelations(
304 B_0 + 0.5 * B_N, G, grad_u, u, u_prev, x_position, t, dt,
305 constitutive_setting, medium, this->current_states_[ip],
306 this->prev_states_[ip], this->material_states_[ip],
307 this->output_data_[ip]);
308
309 auto const B = B_0 + B_N;
310
311 auto const& sigma = this->current_states_[ip].stress_data.sigma;
312 auto const& b = *CD.volumetric_body_force;
313 auto const& C = CD.s_mech_data.stiffness_tensor;
314
315 local_b.noalias() -=
316 (B.transpose() * sigma - N_u_op(N).transpose() * b) * w;
317
318 auto const sigma_geom =
321
322 local_Jac.noalias() += G.transpose() * sigma_geom * G * w;
323
324 local_Jac.noalias() += B.transpose() * C * B * w;
325 }
326 }
327
328 void postTimestepConcrete(Eigen::VectorXd const& local_x,
329 Eigen::VectorXd const& local_x_prev,
330 double const t, double const dt,
331 int const /*process_id*/) override
332 {
333 unsigned const n_integration_points =
335
337 x_position.setElementID(this->element_.getID());
338
340 constitutive_setting;
341
342 auto& medium =
343 *this->process_data_.media_map.getMedium(this->element_.getID());
344
345 for (unsigned ip = 0; ip < n_integration_points; ip++)
346 {
347 x_position.setIntegrationPoint(ip);
348 auto const& N = _ip_data[ip].N;
349 auto const& dNdx = _ip_data[ip].dNdx;
350 auto const x_coord =
351 NumLib::interpolateXCoordinate<ShapeFunction,
353 this->element_, N);
354
355 // For the 2D case the 33-component is needed (and the four entries
356 // of the non-symmetric matrix); In 3d there are nine entries.
357 GradientMatrixType G(DisplacementDim * DisplacementDim +
358 (DisplacementDim == 2 ? 1 : 0),
359 ShapeFunction::NPOINTS * DisplacementDim);
360 Deformation::computeGMatrix<DisplacementDim,
361 ShapeFunction::NPOINTS>(
362 dNdx, G, this->is_axially_symmetric_, N, x_coord);
363
364 GradientVectorType const grad_u = G * local_x;
365
366 // B for Green-Lagrange strain.
368 DisplacementDim, ShapeFunction::NPOINTS,
370 dNdx, N, x_coord, this->is_axially_symmetric_);
371
372 B.noalias() += 0.5 * NonLinearBMatrix::computeBMatrix<
373 DisplacementDim, ShapeFunction::NPOINTS,
375 dNdx, N, grad_u, x_coord,
377
379 B, G, grad_u, local_x, local_x_prev, x_position, t, dt,
380 constitutive_setting, medium, this->current_states_[ip],
381 this->prev_states_[ip], this->material_states_[ip],
382 this->output_data_[ip]);
383
384 this->material_states_[ip].pushBackState();
385 }
386
387 for (unsigned ip = 0; ip < n_integration_points; ip++)
388 {
389 this->prev_states_[ip] = this->current_states_[ip];
390 }
391 }
392
393 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
394 const unsigned integration_point) const override
395 {
396 auto const& N = _secondary_data.N[integration_point];
397
398 // assumes N is stored contiguously in memory
399 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
400 }
401
402private:
403 static constexpr auto localDOF(std::vector<double> const& x)
404 {
405 return NumLib::localDOF<
407 }
408
409private:
410 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
411
413
414 static const int displacement_size =
415 ShapeFunction::NPOINTS * DisplacementDim;
416};
417
418} // namespace LargeDeformation
419} // 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
LargeDeformationLocalAssembler(MeshLib::Element const &e, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, LargeDeformationProcessData< DisplacementDim > &process_data)
typename GMatricesType::GradientVectorType GradientVectorType
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
typename ShapeMatricesType::NodalVectorType NodalVectorType
static constexpr auto localDOF(std::vector< double > const &x)
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
ConstitutiveRelations::ConstitutiveData< DisplacementDim > updateConstitutiveRelations(BMatrixType const &B, GradientMatrixType const &G, GradientVectorType const &grad_u, 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, 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
typename BMatricesType::NodalForceVectorType NodalDisplacementVectorType
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
void postTimestepConcrete(Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, double const t, double const dt, int const) override
typename BMatricesType::NodalForceVectorType NodalForceVectorType
typename GMatricesType::GradientMatrixType GradientMatrixType
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
typename ShapeMatricesType::NodalMatrixType NodalMatrixType
LargeDeformationLocalAssembler(LargeDeformationLocalAssembler &&)=delete
typename BMatricesType::StiffnessMatrixType StiffnessMatrixType
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
LargeDeformationLocalAssembler(LargeDeformationLocalAssembler const &)=delete
constexpr int tensorSize(int dim)
See Tensor type for details.
Definition Tensor.h:19
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)
Eigen::Matrix< double, 3, 3 > kelvinVectorToTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
double determinant(Eigen::MatrixBase< Derived > const &tensor)
Computes determinant of a vectorized tensor.
constexpr Eigen::CwiseNullaryOp< EigenBlockMatrixViewFunctor< D, M >, typename EigenBlockMatrixViewFunctor< D, M >::Matrix > eigenBlockMatrixView(const Eigen::MatrixBase< M > &matrix)
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
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)
void computeGMatrix(DNDX_Type const &dNdx, GMatrixType &g_matrix, const bool is_axially_symmetric, N_Type const &N, const double radius)
Fills a G-matrix based on given shape function dN/dx values.
Definition GMatrix.h:25
Eigen::Matrix< double, MPL::tensorSize(DisplacementDim), MPL::tensorSize(DisplacementDim)> computeSigmaGeom(Eigen::Matrix3d const &sigma_tensor)
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.
BMatrixType computeBMatrix(DNDX_Type const &dNdx, N_Type const &N, GradientVectorType const &grad_u, const double radius, const bool is_axially_symmetric)
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
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, DeformationGradientData< DisplacementDim > const &deformation_gradient_data, GradientVectorType const &deformation_gradient_prev, StatefulData< DisplacementDim > &state, StatefulDataPrev< DisplacementDim > const &prev_state, MaterialStateData< DisplacementDim > &mat_state, ConstitutiveTempData< DisplacementDim > &tmp, ConstitutiveData< DisplacementDim > &cd) const
Evaluate the constitutive setting.
Data that is needed for output purposes, but not directly for the assembly.
DeformationGradientData< DisplacementDim > deformation_gradient_data
Data whose state must be tracked by the process.
std::vector< MaterialStateData< DisplacementDim > > material_states_
std::vector< typename ConstitutiveRelations::StatefulDataPrev< DisplacementDim > > prev_states_
std::vector< typename ConstitutiveRelations::StatefulData< DisplacementDim > > current_states_
std::vector< typename ConstitutiveRelations::OutputData< DisplacementDim > > output_data_
ConstitutiveRelations::SolidConstitutiveRelation< DisplacementDim > const & solid_material_
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N