OGS
SmallDeformationLocalAssemblerMatrix.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <vector>
14
22#include "SecondaryData.h"
24
25namespace ProcessLib
26{
27namespace LIE
28{
29namespace SmallDeformation
30{
31namespace MPL = MaterialPropertyLib;
32
33template <typename ShapeFunction, int DisplacementDim>
36{
37public:
44
49
54
56 MeshLib::Element const& e,
57 std::size_t const local_matrix_size,
58 NumLib::GenericIntegrationMethod const& integration_method,
59 bool const is_axially_symmetric,
61
62 void assemble(double const /*t*/, double const /*dt*/,
63 std::vector<double> const& /*local_x*/,
64 std::vector<double> const& /*local_x_prev*/,
65 std::vector<double>& /*local_M_data*/,
66 std::vector<double>& /*local_K_data*/,
67 std::vector<double>& /*local_b_data*/) override
68 {
70 "SmallDeformationLocalAssembler: assembly without jacobian is not "
71 "implemented.");
72 }
73
74 void assembleWithJacobian(double const t, double const dt,
75 std::vector<double> const& local_x,
76 std::vector<double> const& /*local_x_prev*/,
77 std::vector<double>& local_b_data,
78 std::vector<double>& local_Jac_data) override;
79
80 void preTimestepConcrete(std::vector<double> const& /*local_x*/,
81 double const /*t*/,
82 double const /*delta_t*/) override
83 {
84 unsigned const n_integration_points =
86
87 for (unsigned ip = 0; ip < n_integration_points; ip++)
88 {
89 _ip_data[ip].pushBackState();
90 }
91 }
92
94 double const t, Eigen::VectorXd const& local_x) override;
95
96 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
97 const unsigned integration_point) const override
98 {
99 auto const& N = _secondary_data.N[integration_point];
100
101 // assumes N is stored contiguously in memory
102 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
103 }
104
105 std::vector<double> const& getIntPtSigma(
106 const double t,
107 std::vector<GlobalVector*> const& x,
108 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
109 std::vector<double>& cache) const override;
110
111 std::vector<double> const& getIntPtEpsilon(
112 const double t,
113 std::vector<GlobalVector*> const& x,
114 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
115 std::vector<double>& cache) const override;
116
117 std::vector<double> const& getIntPtFractureStress(
118 const double /*t*/,
119 std::vector<GlobalVector*> const& /*x*/,
120 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
121 std::vector<double>& cache) const override
122 {
123 cache.resize(0);
124 return cache;
125 }
126
127 std::vector<double> const& getIntPtFractureAperture(
128 const double /*t*/,
129 std::vector<GlobalVector*> const& /*x*/,
130 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
131 std::vector<double>& cache) const override
132 {
133 cache.resize(0);
134 return cache;
135 }
136
137private:
139
142 DisplacementDim>;
143 std::vector<IntegrationPointDataType,
144 Eigen::aligned_allocator<IntegrationPointDataType>>
146
151};
152
153} // namespace SmallDeformation
154} // namespace LIE
155} // namespace ProcessLib
156
#define OGS_FATAL(...)
Definition Error.h:26
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
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