OGS
HydroMechanicsLocalAssemblerFracture.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
13
14namespace ProcessLib
15{
16namespace LIE
17{
18namespace HydroMechanics
19{
20template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
21
22 int DisplacementDim>
25{
26public:
31
33 MeshLib::Element const& e,
34 std::size_t const local_matrix_size,
35 std::vector<unsigned> const& dofIndex_to_localIndex,
36 NumLib::GenericIntegrationMethod const& integration_method,
37 bool const is_axially_symmetric,
39
40 void preTimestepConcrete(std::vector<double> const& /*local_x*/,
41 double const /*t*/,
42 double const /*delta_t*/) override
43 {
44 for (auto& data : _ip_data)
45 {
46 data.pushBackState();
47 }
48 }
49
51 const double t, double const dt,
52 Eigen::VectorXd const& local_x) override;
53
54 std::vector<double> const& getIntPtSigma(
55 const double /*t*/,
56 std::vector<GlobalVector*> const& /*x*/,
57 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
58 std::vector<double>& cache) const override
59 {
60 cache.resize(0);
61 return cache;
62 }
63
64 std::vector<double> const& getIntPtEpsilon(
65 const double /*t*/,
66 std::vector<GlobalVector*> const& /*x*/,
67 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
68 std::vector<double>& cache) const override
69 {
70 cache.resize(0);
71 return cache;
72 }
73
74 std::vector<double> const& getIntPtDarcyVelocity(
75 const double /*t*/,
76 std::vector<GlobalVector*> const& /*x*/,
77 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
78 std::vector<double>& cache) const override
79 {
80 cache.resize(0);
81 return cache;
82 }
83
84 std::vector<double> const& getIntPtFractureVelocity(
85 const double t,
86 std::vector<GlobalVector*> const& x,
87 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
88 std::vector<double>& cache) const override;
89
90 std::vector<double> const& getIntPtFractureStress(
91 const double t,
92 std::vector<GlobalVector*> const& x,
93 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
94 std::vector<double>& cache) const override;
95
96 std::vector<double> const& getIntPtFractureAperture(
97 const double t,
98 std::vector<GlobalVector*> const& x,
99 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
100 std::vector<double>& cache) const override;
101
102 std::vector<double> const& getIntPtFracturePermeability(
103 const double t,
104 std::vector<GlobalVector*> const& x,
105 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
106 std::vector<double>& cache) const override;
107
108 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
109 const unsigned integration_point) const override
110 {
111 auto const& N = _secondary_data.N[integration_point];
112
113 // assumes N is stored contiguously in memory
114 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
115 }
116
117private:
118 void assembleWithJacobianConcrete(double const t, double const dt,
119 Eigen::VectorXd const& local_x,
120 Eigen::VectorXd const& local_x_prev,
121 Eigen::VectorXd& local_b,
122 Eigen::MatrixXd& local_J) override;
123
125 double const t, double const dt,
126 Eigen::Ref<const Eigen::VectorXd> const& p,
127 Eigen::Ref<const Eigen::VectorXd> const& p_prev,
128 Eigen::Ref<const Eigen::VectorXd> const& g,
129 Eigen::Ref<const Eigen::VectorXd> const& g_prev,
130 Eigen::Ref<Eigen::VectorXd> rhs_p, Eigen::Ref<Eigen::VectorXd> rhs_g,
131 Eigen::Ref<Eigen::MatrixXd> J_pp, Eigen::Ref<Eigen::MatrixXd> J_pg,
132 Eigen::Ref<Eigen::MatrixXd> J_gg, Eigen::Ref<Eigen::MatrixXd> J_gp);
133
134 // Types for displacement.
140
141 // Types for pressure.
144
145 // Types for the integration point data
148 DisplacementDim>;
149
150 using GlobalDimVectorType = Eigen::Matrix<double, DisplacementDim, 1>;
151
152private:
155
156 std::vector<IntegrationPointDataType,
157 Eigen::aligned_allocator<IntegrationPointDataType>>
159
160 static const int pressure_index = 0;
161 static const int pressure_size = ShapeFunctionPressure::NPOINTS;
162 static const int displacement_index = ShapeFunctionPressure::NPOINTS;
163 static const int displacement_size =
164 ShapeFunctionDisplacement::NPOINTS * DisplacementDim;
165
167 typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType>
169};
170
171} // namespace HydroMechanics
172} // namespace LIE
173} // namespace ProcessLib
174
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
std::vector< double > const & getIntPtFractureVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
HMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > HMatricesType
void postTimestepConcreteWithVector(const double t, double const dt, Eigen::VectorXd const &local_x) override
std::vector< double > const & getIntPtEpsilon(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
void assembleWithJacobianConcrete(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, Eigen::VectorXd &local_b, Eigen::MatrixXd &local_J) override
HydroMechanicsLocalAssemblerFracture(HydroMechanicsLocalAssemblerFracture &&)=delete
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
HydroMechanicsLocalAssemblerFracture(HydroMechanicsLocalAssemblerFracture const &)=delete
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
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
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
void preTimestepConcrete(std::vector< double > const &, double const, double const) override
std::vector< double > const & getIntPtDarcyVelocity(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
IntegrationPointDataFracture< HMatricesType, ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim > IntegrationPointDataType
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< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
void assembleBlockMatricesWithJacobian(double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &p_prev, Eigen::Ref< const Eigen::VectorXd > const &g, Eigen::Ref< const Eigen::VectorXd > const &g_prev, Eigen::Ref< Eigen::VectorXd > rhs_p, Eigen::Ref< Eigen::VectorXd > rhs_g, Eigen::Ref< Eigen::MatrixXd > J_pp, Eigen::Ref< Eigen::MatrixXd > J_pg, Eigen::Ref< Eigen::MatrixXd > J_gg, Eigen::Ref< Eigen::MatrixXd > J_gp)
std::vector< double > const & getIntPtFracturePermeability(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
HydroMechanicsLocalAssemblerInterface(MeshLib::Element const &element, bool const is_axially_symmetric, NumLib::GenericIntegrationMethod const &integration_method, std::size_t n_local_size, std::vector< unsigned > dofIndex_to_localIndex)
Used for the extrapolation of the integration point values.