OGS
SmallDeformationLocalAssemblerMatrix.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <optional>
14#include <vector>
15
23#include "SecondaryData.h"
25
26namespace ProcessLib
27{
28namespace LIE
29{
30namespace SmallDeformation
31{
32namespace MPL = MaterialPropertyLib;
33
34template <typename ShapeFunction, int DisplacementDim>
37{
38public:
46
51
56
58 MeshLib::Element const& e,
59 std::size_t const local_matrix_size,
60 NumLib::GenericIntegrationMethod const& integration_method,
61 bool const is_axially_symmetric,
63
64 void assemble(double const /*t*/, double const /*dt*/,
65 std::vector<double> const& /*local_x*/,
66 std::vector<double> const& /*local_x_prev*/,
67 std::vector<double>& /*local_M_data*/,
68 std::vector<double>& /*local_K_data*/,
69 std::vector<double>& /*local_b_data*/) override
70 {
72 "SmallDeformationLocalAssembler: assembly without jacobian is not "
73 "implemented.");
74 }
75
76 void assembleWithJacobian(double const t, double const dt,
77 std::vector<double> const& local_x,
78 std::vector<double> const& /*local_x_prev*/,
79 std::vector<double>& local_b_data,
80 std::vector<double>& local_Jac_data) override;
81
82 void preTimestepConcrete(std::vector<double> const& /*local_x*/,
83 double const /*t*/,
84 double const /*delta_t*/) override
85 {
86 unsigned const n_integration_points =
88
89 for (unsigned ip = 0; ip < n_integration_points; ip++)
90 {
91 _ip_data[ip].pushBackState();
92 }
93 }
94
96 double const t, Eigen::VectorXd const& local_x) override;
97
98 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
99 const unsigned integration_point) const override
100 {
101 auto const& N = _secondary_data.N[integration_point];
102
103 // assumes N is stored contiguously in memory
104 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
105 }
106
107 std::vector<double> const& getIntPtSigma(
108 const double t,
109 std::vector<GlobalVector*> const& x,
110 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
111 std::vector<double>& cache) const override;
112
113 std::vector<double> const& getIntPtEpsilon(
114 const double t,
115 std::vector<GlobalVector*> const& x,
116 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
117 std::vector<double>& cache) const override;
118
119 std::vector<double> const& getIntPtFractureStress(
120 const double /*t*/,
121 std::vector<GlobalVector*> const& /*x*/,
122 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
123 std::vector<double>& cache) const override
124 {
125 cache.resize(0);
126 return cache;
127 }
128
129 std::vector<double> const& getIntPtFractureAperture(
130 const double /*t*/,
131 std::vector<GlobalVector*> const& /*x*/,
132 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
133 std::vector<double>& cache) const override
134 {
135 cache.resize(0);
136 return cache;
137 }
138
139private:
141
144 DisplacementDim>;
145 std::vector<IntegrationPointDataType,
146 Eigen::aligned_allocator<IntegrationPointDataType>>
148
153
154 std::optional<BBarMatrixType> getDilatationalBBarMatrix() const
155 {
156 if (!(_process_data.use_b_bar))
157 {
158 return std::nullopt;
159 }
160
162 DisplacementDim, ShapeFunction::NPOINTS, ShapeFunction,
166 }
167};
168
169} // namespace SmallDeformation
170} // namespace LIE
171} // namespace ProcessLib
172
#define OGS_FATAL(...)
Definition Error.h:26
MatrixType< 3, ShapeFunction::NPOINTS > BBarMatrixType
MatrixType< _number_of_dof, _number_of_dof > StiffnessMatrixType
VectorType< _number_of_dof > NodalForceVectorType
Rhs residual.
void computeSecondaryVariableConcreteWithVector(double const t, Eigen::VectorXd const &local_x) override
std::vector< double > const & getIntPtFractureStress(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
std::vector< double > const & getIntPtFractureAperture(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
SmallDeformationLocalAssemblerMatrix(SmallDeformationLocalAssemblerMatrix &&)=delete
IntegrationPointDataMatrix< ShapeMatricesType, BMatricesType, DisplacementDim > IntegrationPointDataType
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
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
void preTimestepConcrete(std::vector< double > const &, double const, double const) override
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) 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
BBarMatrixType computeDilatationalBbar(std::vector< IpData, Eigen::aligned_allocator< IpData > > const &ip_data, MeshLib::Element const &element, NumLib::GenericIntegrationMethod const &integration_method, const bool is_axially_symmetric)
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N