OGS
HydroMechanicsLocalAssemblerInterface.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <utility>
14#include <vector>
15
16#include "BaseLib/Error.h"
21
22namespace ProcessLib
23{
24namespace LIE
25{
26namespace HydroMechanics
27{
28
30template <typename ShapeMatrixType>
32{
33 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N;
34};
35
39{
40public:
42 MeshLib::Element const& element,
43 bool const is_axially_symmetric,
44 NumLib::GenericIntegrationMethod const& integration_method,
45 std::size_t n_local_size,
46 std::vector<unsigned>
47 dofIndex_to_localIndex)
48 : _element(element),
49 _is_axially_symmetric(is_axially_symmetric),
50 _integration_method(integration_method),
51 _dofIndex_to_localIndex(std::move(dofIndex_to_localIndex))
52 {
53 _local_u.resize(n_local_size);
54 _local_u_prev.resize(n_local_size);
55 _local_b.resize(_local_u.size());
56 _local_J.resize(_local_u.size(), _local_u.size());
57 }
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_rhs_data*/) override
65 {
67 "HydroMechanicsLocalAssembler: assembly without jacobian is not "
68 "implemented.");
69 }
70
71 void assembleWithJacobian(double const t, double const dt,
72 std::vector<double> const& local_x_,
73 std::vector<double> const& local_x_prev_,
74 std::vector<double>& local_b_data,
75 std::vector<double>& local_Jac_data) override
76 {
77 auto const local_dof_size = local_x_.size();
78
79 _local_u.setZero();
80 for (unsigned i = 0; i < local_dof_size; i++)
81 {
82 _local_u[_dofIndex_to_localIndex[i]] = local_x_[i];
83 }
84 _local_u_prev.setZero();
85 for (unsigned i = 0; i < local_dof_size; i++)
86 {
87 _local_u_prev[_dofIndex_to_localIndex[i]] = local_x_prev_[i];
88 }
89 _local_b.setZero();
90 _local_J.setZero();
91
93 _local_J);
94
95 local_b_data.resize(local_dof_size);
96 for (unsigned i = 0; i < local_dof_size; i++)
97 {
98 local_b_data[i] = _local_b[_dofIndex_to_localIndex[i]];
99 }
100
101 local_Jac_data.resize(local_dof_size * local_dof_size);
102 for (unsigned i = 0; i < local_dof_size; i++)
103 {
104 for (unsigned j = 0; j < local_dof_size; j++)
105 {
106 local_Jac_data[i * local_dof_size + j] = _local_J(
108 }
109 }
110 }
111
112 void postTimestepConcrete(Eigen::VectorXd const& local_x_,
113 Eigen::VectorXd const& /*local_x_prev*/,
114 const double t, double const dt,
115 int const /*process_id*/) override
116 {
117 _local_u.setZero();
118 for (Eigen::Index i = 0; i < local_x_.rows(); i++)
119 {
120 _local_u[_dofIndex_to_localIndex[i]] = local_x_[i];
121 }
122
124 }
125
126 virtual std::vector<double> const& getIntPtSigma(
127 const double t,
128 std::vector<GlobalVector*> const& x,
129 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
130 std::vector<double>& cache) const = 0;
131
132 virtual std::vector<double> const& getIntPtEpsilon(
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 = 0;
137
138 virtual std::vector<double> const& getIntPtDarcyVelocity(
139 const double t,
140 std::vector<GlobalVector*> const& x,
141 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
142 std::vector<double>& cache) const = 0;
143
144 virtual std::vector<double> const& getIntPtFractureVelocity(
145 const double t,
146 std::vector<GlobalVector*> const& x,
147 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
148 std::vector<double>& cache) const = 0;
149
150 virtual std::vector<double> const& getIntPtFractureStress(
151 const double t,
152 std::vector<GlobalVector*> const& x,
153 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
154 std::vector<double>& cache) const = 0;
155
156 virtual std::vector<double> const& getIntPtFractureAperture(
157 const double t,
158 std::vector<GlobalVector*> const& x,
159 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
160 std::vector<double>& cache) const = 0;
161
162 virtual std::vector<double> const& getIntPtFracturePermeability(
163 const double t,
164 std::vector<GlobalVector*> const& x,
165 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
166 std::vector<double>& cache) const = 0;
167
168protected:
170 double const t, double const dt, Eigen::VectorXd const& local_u,
171 Eigen::VectorXd const& local_u_prev, Eigen::VectorXd& local_b,
172 Eigen::MatrixXd& local_J) = 0;
173
175 double const t, double const dt, Eigen::VectorXd const& local_u) = 0;
176
180
181private:
182 Eigen::VectorXd _local_u;
183 Eigen::VectorXd _local_u_prev;
184 Eigen::VectorXd _local_b;
185 Eigen::MatrixXd _local_J;
186 // a vector for mapping the index in the local DoF vector to the index in
187 // the complete local solution vector which also include nodes where DoF are
188 // not assigned.
189 std::vector<unsigned> const _dofIndex_to_localIndex;
190};
191
192} // namespace HydroMechanics
193} // namespace LIE
194} // namespace ProcessLib
Definition of the Element class.
#define OGS_FATAL(...)
Definition Error.h:26
virtual void postTimestepConcreteWithVector(double const t, double const dt, Eigen::VectorXd const &local_u)=0
virtual 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 =0
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)
virtual void assembleWithJacobianConcrete(double const t, double const dt, Eigen::VectorXd const &local_u, Eigen::VectorXd const &local_u_prev, Eigen::VectorXd &local_b, Eigen::MatrixXd &local_J)=0
virtual std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual 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 =0
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x_, std::vector< double > const &local_x_prev_, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
virtual 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 =0
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
virtual 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 =0
void postTimestepConcrete(Eigen::VectorXd const &local_x_, Eigen::VectorXd const &, const double t, double const dt, int const) override
virtual 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 =0
virtual 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 =0
Used for the extrapolation of the integration point values.
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N