OGS
SmallDeformationLocalAssemblerMatrixNearFracture.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <unordered_map>
14#include <vector>
15
25#include "SecondaryData.h"
27
28namespace ProcessLib
29{
30namespace LIE
31{
32namespace SmallDeformation
33{
34namespace MPL = MaterialPropertyLib;
35
36template <typename ShapeFunction, int DisplacementDim>
39{
40public:
48
54
59
61 MeshLib::Element const& e,
62 std::size_t const n_variables,
63 std::size_t const local_matrix_size,
64 std::vector<unsigned> const& dofIndex_to_localIndex,
65 NumLib::GenericIntegrationMethod const& integration_method,
66 bool const is_axially_symmetric,
68
69 void assemble(double const /*t*/, double const /*dt*/,
70 std::vector<double> const& /*local_x*/,
71 std::vector<double> const& /*local_x_prev*/,
72 std::vector<double>& /*local_M_data*/,
73 std::vector<double>& /*local_K_data*/,
74 std::vector<double>& /*local_b_data*/) override
75 {
77 "SmallDeformationLocalAssembler: assembly without jacobian is not "
78 "implemented.");
79 }
80
81 void assembleWithJacobian(double const t, double const dt,
82 Eigen::VectorXd const& local_u,
83 Eigen::VectorXd& local_b,
84 Eigen::MatrixXd& local_J) override;
85
86 void preTimestepConcrete(std::vector<double> const& /*local_x*/,
87 double const /*t*/,
88 double const /*delta_t*/) override
89 {
90 unsigned const n_integration_points =
92
93 for (unsigned ip = 0; ip < n_integration_points; ip++)
94 {
95 _ip_data[ip].pushBackState();
96 }
97 }
98
100 double const t, Eigen::VectorXd const& local_x) override;
101
102 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
103 const unsigned integration_point) const override
104 {
105 auto const& N = _secondary_data.N[integration_point];
106
107 // assumes N is stored contiguously in memory
108 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
109 }
110
111 std::vector<double> const& getIntPtSigma(
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& getIntPtEpsilon(
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 std::vector<double> const& getIntPtFractureStress(
124 const double /*t*/,
125 std::vector<GlobalVector*> const& /*x*/,
126 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
127 std::vector<double>& cache) const override
128 {
129 cache.resize(0);
130 return cache;
131 }
132
133 std::vector<double> const& getIntPtFractureAperture(
134 const double /*t*/,
135 std::vector<GlobalVector*> const& /*x*/,
136 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
137 std::vector<double>& cache) const override
138 {
139 cache.resize(0);
140 return cache;
141 }
142
143private:
145 std::vector<FractureProperty*> _fracture_props;
146 std::vector<JunctionProperty*> _junction_props;
147 std::unordered_map<int, int> _fracID_to_local;
148
151 DisplacementDim>;
152
153 std::vector<IntegrationPointDataType,
154 Eigen::aligned_allocator<IntegrationPointDataType>>
156
158
162
163 std::optional<BBarMatrixType> getDilatationalBBarMatrix() const
164 {
165 if (!(_process_data.use_b_bar))
166 {
167 return std::nullopt;
168 }
169
171 DisplacementDim, ShapeFunction::NPOINTS, ShapeFunction,
175 }
176};
177
178} // namespace SmallDeformation
179} // namespace LIE
180} // namespace ProcessLib
181
#define OGS_FATAL(...)
Definition Error.h:26
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
MatrixType< 3, ShapeFunction::NPOINTS > BBarMatrixType
MatrixType< _number_of_dof, _number_of_dof > StiffnessMatrixType
VectorType< _number_of_dof > NodalForceVectorType
Rhs residual.
IntegrationPointDataMatrix< ShapeMatricesType, BMatricesType, DisplacementDim > IntegrationPointDataType
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
SmallDeformationLocalAssemblerMatrixNearFracture(SmallDeformationLocalAssemblerMatrixNearFracture const &)=delete
void assembleWithJacobian(double const t, double const dt, Eigen::VectorXd const &local_u, Eigen::VectorXd &local_b, Eigen::MatrixXd &local_J) override
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
SmallDeformationLocalAssemblerMatrixNearFracture(SmallDeformationLocalAssemblerMatrixNearFracture &&)=delete
std::vector< double > const & getIntPtFractureStress(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) 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< 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
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
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