OGS
LargeDeformationFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <limits>
14#include <memory>
15#include <optional>
16#include <vector>
17
24#include "NumLib/DOF/LocalDOF.h"
38
39namespace ProcessLib
40{
41namespace LargeDeformation
42{
43namespace MPL = MaterialPropertyLib;
44
47template <typename ShapeMatrixType>
49{
50 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N;
51};
52
53template <int DisplacementDim, typename ShapeMatricesType>
54Eigen::Matrix<double, MPL::tensorSize(DisplacementDim),
55 MPL::tensorSize(DisplacementDim)>
56computeSigmaGeom(Eigen::Matrix3d const& sigma_tensor)
57{
58 static constexpr auto& sigma_geom_op = MathLib::eigenBlockMatrixView<
59 DisplacementDim,
60 Eigen::Matrix<double, DisplacementDim, DisplacementDim>>;
61
62 using SigmaGeom = Eigen::Matrix<double, MPL::tensorSize(DisplacementDim),
63 MPL::tensorSize(DisplacementDim)>;
64 if constexpr (DisplacementDim == 2)
65 {
66 SigmaGeom sigma_geom = SigmaGeom::Zero(5, 5);
67 sigma_geom.template block<4, 4>(0, 0) =
68 sigma_geom_op(sigma_tensor.template block<2, 2>(0, 0).eval());
69 sigma_geom(4, 4) = sigma_tensor(2, 2);
70
71 return sigma_geom;
72 }
73
74 if constexpr (DisplacementDim == 3)
75 {
76 return sigma_geom_op(sigma_tensor);
77 }
78}
79
80template <typename ShapeFunction, int DisplacementDim>
82 : public LargeDeformationLocalAssemblerInterface<DisplacementDim>
83{
84public:
93
99
103
105
106 using IpData =
108
109 static constexpr auto& N_u_op = MathLib::eigenBlockMatrixView<
110 DisplacementDim, typename ShapeMatricesType::NodalRowVectorType>;
111
113 delete;
115
117 MeshLib::Element const& e,
118 std::size_t const /*local_matrix_size*/,
119 NumLib::GenericIntegrationMethod const& integration_method,
120 bool const is_axially_symmetric,
123 e, integration_method, is_axially_symmetric, process_data)
124 {
125 unsigned const n_integration_points =
127
128 _ip_data.resize(n_integration_points);
129 _secondary_data.N.resize(n_integration_points);
130
131 auto const shape_matrices =
133 DisplacementDim>(
134 e, is_axially_symmetric, this->integration_method_);
135
136 for (unsigned ip = 0; ip < n_integration_points; ip++)
137 {
138 auto& ip_data = _ip_data[ip];
139 auto const& sm = shape_matrices[ip];
140 _ip_data[ip].integration_weight =
142 sm.integralMeasure * sm.detJ;
143
144 ip_data.N_u = sm.N;
145 ip_data.dNdx_u = sm.dNdx;
146
147 _secondary_data.N[ip] = shape_matrices[ip].N;
148 }
149 }
150
151 void initializeConcrete() override
152 {
153 unsigned const n_integration_points =
155 for (unsigned ip = 0; ip < n_integration_points; ip++)
156 {
157 auto const& ip_data = _ip_data[ip];
158
159 ParameterLib::SpatialPosition const x_position{
160 std::nullopt, this->element_.getID(),
162 NumLib::interpolateCoordinates<ShapeFunction,
164 this->element_, ip_data.N_u))};
165
167 if (this->process_data_.initial_stress != nullptr)
168 {
169 this->current_states_[ip].stress_data.sigma.noalias() =
171 DisplacementDim>((*this->process_data_.initial_stress)(
172 std::numeric_limits<
173 double>::quiet_NaN() /* time independent */,
174 x_position));
175 }
176
177 double const t = 0; // TODO (naumov) pass t from top
178 auto& material_state = this->material_states_[ip];
179 this->solid_material_.initializeInternalStateVariables(
180 t, x_position, *material_state.material_state_variables);
181
182 material_state.pushBackState();
183 this->prev_states_[ip] = this->current_states_[ip];
184 }
185 }
186
188 double const detF0,
189 BMatrixType const&
190 B /*B_{linear}+B_{nonlinear}/2 for Green-Lagrange strain*/,
191 GradientVectorType const& grad_u,
192 Eigen::Ref<Eigen::VectorXd const> const& u,
194 output_data) const
195 {
196 // Note: Here B=B_{linear}+B_{nonlinear}/2 For Green-Lagrange strain.
197 output_data.eps_data.eps = B * u;
198 output_data.deformation_gradient_data.deformation_gradient =
200 output_data.deformation_gradient_data.volume_ratio =
202 output_data.deformation_gradient_data.deformation_gradient);
203
204 if (this->process_data_.bar_det_f_type ==
206 {
207 return 1.0;
208 }
209
210 double const detF = output_data.deformation_gradient_data.volume_ratio;
211 double const J_ratio = detF0 / detF;
212
213 if (J_ratio < .0)
214 {
215 OGS_FATAL(
216 "det(F0)/det(F) is negative with det(F0) = {:g}, and det(F) = "
217 "{:g} ",
218 detF0, detF);
219 }
220
221 double const alpha =
222 DisplacementDim == 3 ? std::cbrt(J_ratio) : std::sqrt(J_ratio);
223
224 output_data.deformation_gradient_data.deformation_gradient
225 .template segment<DisplacementDim * DisplacementDim>(0) *= alpha;
226 output_data.deformation_gradient_data.volume_ratio *=
227 std::pow(alpha, DisplacementDim);
228
229 double const alpha_p2 = alpha * alpha;
230 output_data.eps_data.eps =
231 alpha_p2 * output_data.eps_data.eps +
232 0.5 * (alpha_p2 - 1) *
235 return alpha;
236 }
237
240 double const alpha, GradientMatrixType const& G,
241 Eigen::Ref<Eigen::VectorXd const> const& u_prev,
242 ParameterLib::SpatialPosition const& x_position, double const t,
243 double const dt,
245 CS,
246 MaterialPropertyLib::Medium const& medium,
248 current_state,
250 prev_state,
253 output_data) const
254 {
255 double const T_ref =
256 this->process_data_.reference_temperature
257 ? (*this->process_data_.reference_temperature)(t, x_position)[0]
258 : std::numeric_limits<double>::quiet_NaN();
259
261 models(this->process_data_, this->solid_material_);
263 tmp;
265
266 CS.eval(
267 models, t, dt, x_position, //
268 medium, //
269 T_ref, //
270 output_data.deformation_gradient_data, //
271 alpha * (G * u_prev +
273 current_state, prev_state, material_state, tmp, CD);
274
275 return CD;
276 }
277
278 void assemble(double const /*t*/, double const /*dt*/,
279 std::vector<double> const& /*local_x*/,
280 std::vector<double> const& /*local_x_prev*/,
281 std::vector<double>& /*local_M_data*/,
282 std::vector<double>& /*local_K_data*/,
283 std::vector<double>& /*local_b_data*/) override
284 {
285 OGS_FATAL(
286 "LargeDeformationLocalAssembler: assembly without jacobian is not "
287 "implemented.");
288 }
289
290 std::optional<std::tuple<double, VectorTypeForFbar>> computeFBarVariables(
291 bool const compute_detF0_only, Eigen::VectorXd const& local_x) const
292 {
293 if (this->process_data_.bar_det_f_type ==
295 {
296 return std::nullopt;
297 }
298
299 if (this->process_data_.bar_det_f_type ==
301 {
303 DisplacementDim, GradientVectorType, VectorTypeForFbar,
304 NodalVectorType, ShapeFunction, ShapeMatricesType, IpData>(
305 _ip_data, compute_detF0_only, local_x,
306 this->integration_method_, this->element_,
308 }
309
311 DisplacementDim, GradientVectorType, GradientMatrixType,
312 VectorTypeForFbar, ShapeFunction, ShapeMatricesType>(
313 compute_detF0_only, local_x, this->element_,
315 }
316
317 void assembleWithJacobian(double const t, double const dt,
318 std::vector<double> const& local_x,
319 std::vector<double> const& local_x_prev,
320 std::vector<double>& local_b_data,
321 std::vector<double>& local_Jac_data) override
322 {
323 auto const local_matrix_size = local_x.size();
324
326 local_Jac_data, local_matrix_size, local_matrix_size);
327
329 local_b_data, local_matrix_size);
330
331 auto [u] = localDOF(local_x);
332 auto [u_prev] = localDOF(local_x_prev);
333
334 unsigned const n_integration_points =
336
338 x_position.setElementID(this->element_.getID());
339
341 constitutive_setting;
342 auto const& medium =
343 *this->process_data_.media_map.getMedium(this->element_.getID());
344
345 auto const f_bar_variables =
346 computeFBarVariables(/*compute_detF0_only?*/ false, u);
347
348 for (unsigned ip = 0; ip < n_integration_points; ip++)
349 {
350 auto const& w = _ip_data[ip].integration_weight;
351 auto const& N = _ip_data[ip].N_u;
352 auto const& dNdx = _ip_data[ip].dNdx_u;
353
354 auto const x_coord =
355 NumLib::interpolateXCoordinate<ShapeFunction,
357 this->element_, N);
358
359 // For the 2D case the 33-component is needed (and the four entries
360 // of the non-symmetric matrix); In 3d there are nine entries.
361 GradientMatrixType G(DisplacementDim * DisplacementDim +
362 (DisplacementDim == 2 ? 1 : 0),
363 ShapeFunction::NPOINTS * DisplacementDim);
364 Deformation::computeGMatrix<DisplacementDim,
365 ShapeFunction::NPOINTS>(
366 dNdx, G, this->is_axially_symmetric_, N, x_coord);
367
368 GradientVectorType const grad_u = G * u;
369
370 auto const B_0 = LinearBMatrix::computeBMatrix<
371 DisplacementDim, ShapeFunction::NPOINTS,
373 dNdx, N, x_coord, this->is_axially_symmetric_);
374
375 auto const B_N = NonLinearBMatrix::computeBMatrix<
376 DisplacementDim, ShapeFunction::NPOINTS,
378 dNdx, N, grad_u, x_coord, this->is_axially_symmetric_);
379
380 double const detF0 =
381 f_bar_variables ? std::get<double>(*f_bar_variables) : 0.0;
382 double const alpha = computeOutputStrainData(
383 detF0, B_0 + 0.5 * B_N, grad_u, u, this->output_data_[ip]);
384
385 auto const CD = updateConstitutiveRelations(
386 alpha, G, u_prev, x_position, t, dt, constitutive_setting,
387 medium, this->current_states_[ip], this->prev_states_[ip],
388 this->material_states_[ip], this->output_data_[ip]);
389
390 BMatrixType B = B_0 + B_N;
391
392 auto const& sigma = this->current_states_[ip].stress_data.sigma;
393 auto const& b = *CD.volumetric_body_force;
394 auto const& C = CD.s_mech_data.stiffness_tensor;
395
396 auto const sigma_geom =
399
400 if (!f_bar_variables)
401 {
402 local_b.noalias() -=
403 (B.transpose() * sigma - N_u_op(N).transpose() * b) * w;
404
405 local_Jac.noalias() += G.transpose() * sigma_geom * G * w;
406
407 local_Jac.noalias() += B.transpose() * C * B * w;
408 }
409 else // with F bar
410 {
411 double const alpha_p2 = alpha * alpha;
412 local_Jac.noalias() +=
413 alpha_p2 * G.transpose() * sigma_geom * G * w;
414
415 auto const q0 = std::get<VectorTypeForFbar>(*f_bar_variables);
416
417 auto const& F =
418 this->output_data_[ip]
419 .deformation_gradient_data.deformation_gradient;
421 DisplacementDim, ShapeFunction::NPOINTS, VectorTypeForFbar,
424 dNdx, F, this->is_axially_symmetric_);
425
426 auto const q0_q = q0 - q;
427
428 auto const S_q = sigma * (q0_q).transpose();
429
430 local_Jac.noalias() +=
431 2.0 * alpha_p2 * S_q.transpose() * B * w / DisplacementDim;
432
433 auto const twoEplsI =
435 alpha_p2, this->output_data_[ip].eps_data.eps,
436 this->is_axially_symmetric_);
437
438 // B bar:
439 BMatrixType const Bbar =
440 alpha_p2 *
441 (B + twoEplsI * (q0_q).transpose() / DisplacementDim);
442
443 local_Jac.noalias() +=
444 2.0 * Bbar.transpose() * S_q * w / DisplacementDim;
445
446 local_Jac.noalias() += Bbar.transpose() * C * Bbar * w;
447
448 local_Jac.noalias() -=
449 alpha_p2 * sigma.dot(twoEplsI) *
450 (q0 * q0.transpose() - q * q.transpose()) * w /
451 DisplacementDim;
452
453 local_b.noalias() -=
454 (Bbar.transpose() * sigma - N_u_op(N).transpose() * b) * w;
455 }
456 }
457 }
458
459 void postTimestepConcrete(Eigen::VectorXd const& local_x,
460 Eigen::VectorXd const& local_x_prev,
461 double const t, double const dt,
462 int const /*process_id*/) override
463 {
464 unsigned const n_integration_points =
466
468 x_position.setElementID(this->element_.getID());
469
471 constitutive_setting;
472
473 auto& medium =
474 *this->process_data_.media_map.getMedium(this->element_.getID());
475
476 auto const f_bar_variables =
477 computeFBarVariables(/*compute_detF0_only?*/ true, local_x);
478
479 for (unsigned ip = 0; ip < n_integration_points; ip++)
480 {
481 auto const& N = _ip_data[ip].N_u;
482 auto const& dNdx = _ip_data[ip].dNdx_u;
483 auto const x_coord =
484 NumLib::interpolateXCoordinate<ShapeFunction,
486 this->element_, N);
487
488 // For the 2D case the 33-component is needed (and the four entries
489 // of the non-symmetric matrix); In 3d there are nine entries.
490 GradientMatrixType G(DisplacementDim * DisplacementDim +
491 (DisplacementDim == 2 ? 1 : 0),
492 ShapeFunction::NPOINTS * DisplacementDim);
493 Deformation::computeGMatrix<DisplacementDim,
494 ShapeFunction::NPOINTS>(
495 dNdx, G, this->is_axially_symmetric_, N, x_coord);
496
497 GradientVectorType const grad_u = G * local_x;
498
499 // B for Green-Lagrange strain.
501 DisplacementDim, ShapeFunction::NPOINTS,
503 dNdx, N, x_coord, this->is_axially_symmetric_);
504
505 B.noalias() += 0.5 * NonLinearBMatrix::computeBMatrix<
506 DisplacementDim, ShapeFunction::NPOINTS,
508 dNdx, N, grad_u, x_coord,
510
511 double const detF0 =
512 f_bar_variables ? std::get<double>(*f_bar_variables) : 0.0;
513 double const alpha = computeOutputStrainData(
514 detF0, B, grad_u, local_x, this->output_data_[ip]);
515
517 alpha, G, local_x_prev, x_position, t, dt, constitutive_setting,
518 medium, this->current_states_[ip], this->prev_states_[ip],
519 this->material_states_[ip], this->output_data_[ip]);
520
521 this->material_states_[ip].pushBackState();
522 }
523
524 for (unsigned ip = 0; ip < n_integration_points; ip++)
525 {
526 this->prev_states_[ip] = this->current_states_[ip];
527 }
528 }
529
530 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
531 const unsigned integration_point) const override
532 {
533 auto const& N = _secondary_data.N[integration_point];
534
535 // assumes N is stored contiguously in memory
536 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
537 }
538
539private:
540 static constexpr auto localDOF(std::vector<double> const& x)
541 {
542 return NumLib::localDOF<
544 }
545
546private:
547 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
548
550
551 static const int displacement_size =
552 ShapeFunction::NPOINTS * DisplacementDim;
553};
554
555} // namespace LargeDeformation
556} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
constexpr 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)
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
std::optional< std::tuple< double, VectorTypeForFbar > > computeFBarVariables(bool const compute_detF0_only, Eigen::VectorXd const &local_x) const
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
double computeOutputStrainData(double const detF0, BMatrixType const &B, GradientVectorType const &grad_u, Eigen::Ref< Eigen::VectorXd const > const &u, typename ConstitutiveRelations::OutputData< DisplacementDim > &output_data) const
typename ShapeMatricesType::NodalVectorType NodalVectorType
static constexpr auto localDOF(std::vector< double > const &x)
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
typename BMatricesType::NodalForceVectorType NodalDisplacementVectorType
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
typename BMatricesType::NodalForceVectorType VectorTypeForFbar
void postTimestepConcrete(Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, double const t, double const dt, int const) override
ConstitutiveRelations::ConstitutiveData< DisplacementDim > updateConstitutiveRelations(double const alpha, GradientMatrixType const &G, 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 NodalForceVectorType
typename GMatricesType::GradientMatrixType GradientMatrixType
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
typename ShapeMatricesType::GlobalDimNodalMatrixType GlobalDimNodalMatrixType
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)
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > identityForF(bool const is_axially_symmetric)
std::tuple< double, VectorTypeForFbar > computeFBarInitialVariablesElementCenter(bool const compute_detF0_only, Eigen::VectorXd const &u, MeshLib::Element const &element, bool const is_axially_symmetric)
VectorTypeForFbar computeQVector(DNDX_Type const &dNdx, GradientVectorType const &F, bool const)
std::tuple< double, VectorTypeForFbar > computeFBarInitialVariablesAverage(std::vector< IpData, Eigen::aligned_allocator< IpData > > const &ip_data, bool const compute_detF0_only, Eigen::VectorXd const &u, NumLib::GenericIntegrationMethod const &integration_method, MeshLib::Element const &element, bool const is_axially_symmetric)
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > compute2EPlusI(double const alpha_p2, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps_bar, bool const is_axially_symmetric)
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType
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