OGS
SmallDeformationLocalAssemblerFracture.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{
31template <typename ShapeFunction, int DisplacementDim>
34{
35public:
42
48
50 using GlobalDimVectorType = Eigen::Matrix<double, DisplacementDim, 1>;
51
56
58 MeshLib::Element const& e,
59 std::size_t const n_variables,
60 std::size_t const local_matrix_size,
61 std::vector<unsigned> const& dofIndex_to_localIndex,
62 NumLib::GenericIntegrationMethod const& integration_method,
63 bool const is_axially_symmetric,
65
66 void assemble(double const /*t*/, double const /*dt*/,
67 std::vector<double> const& /*local_x*/,
68 std::vector<double> const& /*local_x_prev*/,
69 std::vector<double>& /*local_M_data*/,
70 std::vector<double>& /*local_K_data*/,
71 std::vector<double>& /*local_b_data*/) override
72 {
74 "SmallDeformationLocalAssembler: assembly without jacobian is not "
75 "implemented.");
76 }
77
78 void assembleWithJacobian(double const t, double const dt,
79 Eigen::VectorXd const& local_u,
80 Eigen::VectorXd& local_b,
81 Eigen::MatrixXd& local_J) override;
82
83 void preTimestepConcrete(std::vector<double> const& /*local_x*/,
84 double const /*t*/,
85 double const /*delta_t*/) override
86 {
87 unsigned const n_integration_points =
89
90 for (unsigned ip = 0; ip < n_integration_points; ip++)
91 {
92 _ip_data[ip].pushBackState();
93 }
94 }
95
97 const double t, Eigen::VectorXd const& local_u) override;
98
99 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
100 const unsigned integration_point) const override
101 {
102 auto const& N = _secondary_data.N[integration_point];
103
104 // assumes N is stored contiguously in memory
105 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
106 }
107
108 std::vector<double> const& getIntPtSigma(
109 const double /*t*/,
110 std::vector<GlobalVector*> const& /*x*/,
111 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
112 std::vector<double>& cache) const override
113 {
114 cache.resize(0);
115 return cache;
116 }
117
118 std::vector<double> const& getIntPtEpsilon(
119 const double /*t*/,
120 std::vector<GlobalVector*> const& /*x*/,
121 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
122 std::vector<double>& cache) const override
123 {
124 cache.resize(0);
125 return cache;
126 }
127
128 std::vector<double> const& getIntPtFractureStress(
129 const double t,
130 std::vector<GlobalVector*> const& x,
131 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
132 std::vector<double>& cache) const override;
133
134 std::vector<double> const& getIntPtFractureAperture(
135 const double t,
136 std::vector<GlobalVector*> const& x,
137 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
138 std::vector<double>& cache) const override;
139
140private:
142 std::vector<FractureProperty*> _fracture_props;
143 std::vector<JunctionProperty*> _junction_props;
144 std::unordered_map<int, int> _fracID_to_local;
146
149 std::vector<IntegrationPointDataType,
150 Eigen::aligned_allocator<IntegrationPointDataType>>
152
154 std::vector<ShapeMatrices, Eigen::aligned_allocator<
159};
160
161} // namespace SmallDeformation
162} // namespace LIE
163} // namespace ProcessLib
164
#define OGS_FATAL(...)
Definition Error.h:26
MatrixType< _number_of_dof, _number_of_dof > StiffnessMatrixType
VectorType< _number_of_dof > NodalForceVectorType
MatrixType< DisplacementDim, _number_of_dof > HMatrixType
VectorType< DisplacementDim > ForceVectorType
IntegrationPointDataFracture< HMatricesType, DisplacementDim > IntegrationPointDataType
std::vector< double > const & getIntPtFractureAperture(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
void computeSecondaryVariableConcreteWithVector(const double t, Eigen::VectorXd const &local_u) override
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
void assembleWithJacobian(double const t, double const dt, Eigen::VectorXd const &local_u, Eigen::VectorXd &local_b, Eigen::MatrixXd &local_J) override
std::vector< double > const & getIntPtEpsilon(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
SmallDeformationLocalAssemblerFracture(SmallDeformationLocalAssemblerFracture &&)=delete
void preTimestepConcrete(std::vector< double > const &, double const, double const) override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
SmallDeformationLocalAssemblerFracture(SmallDeformationLocalAssemblerFracture const &)=delete
std::vector< double > const & getIntPtSigma(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
std::vector< ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > _shape_matrices
std::vector< double > const & getIntPtFractureStress(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
Coordinates mapping matrices at particular location.
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N