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:
47
53
58
60 MeshLib::Element const& e,
61 std::size_t const n_variables,
62 std::size_t const local_matrix_size,
63 std::vector<unsigned> const& dofIndex_to_localIndex,
64 NumLib::GenericIntegrationMethod const& integration_method,
65 bool const is_axially_symmetric,
67
68 void assemble(double const /*t*/, double const /*dt*/,
69 std::vector<double> const& /*local_x*/,
70 std::vector<double> const& /*local_x_prev*/,
71 std::vector<double>& /*local_M_data*/,
72 std::vector<double>& /*local_K_data*/,
73 std::vector<double>& /*local_b_data*/) override
74 {
76 "SmallDeformationLocalAssembler: assembly without jacobian is not "
77 "implemented.");
78 }
79
80 void assembleWithJacobian(double const t, double const dt,
81 Eigen::VectorXd const& local_u,
82 Eigen::VectorXd& local_b,
83 Eigen::MatrixXd& local_J) override;
84
85 void preTimestepConcrete(std::vector<double> const& /*local_x*/,
86 double const /*t*/,
87 double const /*delta_t*/) override
88 {
89 unsigned const n_integration_points =
91
92 for (unsigned ip = 0; ip < n_integration_points; ip++)
93 {
94 _ip_data[ip].pushBackState();
95 }
96 }
97
99 double const t, Eigen::VectorXd const& local_x) override;
100
101 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
102 const unsigned integration_point) const override
103 {
104 auto const& N = _secondary_data.N[integration_point];
105
106 // assumes N is stored contiguously in memory
107 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
108 }
109
110 std::vector<double> const& getIntPtSigma(
111 const double t,
112 std::vector<GlobalVector*> const& x,
113 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
114 std::vector<double>& cache) const override;
115
116 std::vector<double> const& getIntPtEpsilon(
117 const double t,
118 std::vector<GlobalVector*> const& x,
119 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
120 std::vector<double>& cache) const override;
121
122 std::vector<double> const& getIntPtFractureStress(
123 const double /*t*/,
124 std::vector<GlobalVector*> const& /*x*/,
125 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
126 std::vector<double>& cache) const override
127 {
128 cache.resize(0);
129 return cache;
130 }
131
132 std::vector<double> const& getIntPtFractureAperture(
133 const double /*t*/,
134 std::vector<GlobalVector*> const& /*x*/,
135 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
136 std::vector<double>& cache) const override
137 {
138 cache.resize(0);
139 return cache;
140 }
141
142private:
144 std::vector<FractureProperty*> _fracture_props;
145 std::vector<JunctionProperty*> _junction_props;
146 std::unordered_map<int, int> _fracID_to_local;
147
150 DisplacementDim>;
151
152 std::vector<IntegrationPointDataType,
153 Eigen::aligned_allocator<IntegrationPointDataType>>
155
157
161};
162
163} // namespace SmallDeformation
164} // namespace LIE
165} // namespace ProcessLib
166
#define OGS_FATAL(...)
Definition Error.h:26
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
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
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