OGS
SmallDeformationLocalAssemblerFracture.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <vector>
7
15#include "SecondaryData.h"
17
18namespace ProcessLib
19{
20namespace LIE
21{
22namespace SmallDeformation
23{
24template <typename ShapeFunction, int DisplacementDim>
27{
28public:
35
41
43 using GlobalDimVectorType = Eigen::Matrix<double, DisplacementDim, 1>;
44
49
51 MeshLib::Element const& e,
52 std::size_t const n_variables,
53 std::size_t const local_matrix_size,
54 std::vector<unsigned> const& dofIndex_to_localIndex,
55 NumLib::GenericIntegrationMethod const& integration_method,
56 bool const is_axially_symmetric,
58
59 void assemble(double const /*t*/, double const /*dt*/,
60 std::vector<double> const& /*local_x*/,
61 std::vector<double> const& /*local_x_prev*/,
62 std::vector<double>& /*local_M_data*/,
63 std::vector<double>& /*local_K_data*/,
64 std::vector<double>& /*local_b_data*/) override
65 {
67 "SmallDeformationLocalAssembler: assembly without jacobian is not "
68 "implemented.");
69 }
70
71 void assembleWithJacobian(double const t, double const dt,
72 Eigen::VectorXd const& local_u,
73 Eigen::VectorXd& local_b,
74 Eigen::MatrixXd& local_J) override;
75
76 void preTimestepConcrete(std::vector<double> const& /*local_x*/,
77 double const /*t*/,
78 double const /*delta_t*/) override
79 {
80 unsigned const n_integration_points =
81 _integration_method.getNumberOfPoints();
82
83 for (unsigned ip = 0; ip < n_integration_points; ip++)
84 {
85 _ip_data[ip].pushBackState();
86 }
87 }
88
90 const double t, Eigen::VectorXd const& local_u) override;
91
92 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
93 const unsigned integration_point) const override
94 {
95 auto const& N = _secondary_data.N[integration_point];
96
97 // assumes N is stored contiguously in memory
98 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
99 }
100
101 std::vector<double> const& getIntPtSigma(
102 const double /*t*/,
103 std::vector<GlobalVector*> const& /*x*/,
104 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
105 std::vector<double>& cache) const override
106 {
107 cache.resize(0);
108 return cache;
109 }
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 cache.resize(0);
118 return cache;
119 }
120
121 std::vector<double> const& getIntPtFractureStress(
122 const double t,
123 std::vector<GlobalVector*> const& x,
124 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
125 std::vector<double>& cache) const override;
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
133private:
135 std::vector<FractureProperty*> _fracture_props;
136 std::vector<JunctionProperty*> _junction_props;
137 std::unordered_map<int, int> _fracID_to_local;
139
142 std::vector<IntegrationPointDataType,
143 Eigen::aligned_allocator<IntegrationPointDataType>>
145
147 std::vector<ShapeMatrices, Eigen::aligned_allocator<
152};
153
154} // namespace SmallDeformation
155} // namespace LIE
156} // namespace ProcessLib
157
#define OGS_FATAL(...)
Definition Error.h:19
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
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