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"
20
21namespace ProcessLib
22{
23namespace LIE
24{
25namespace HydroMechanics
26{
27
29template <typename ShapeMatrixType>
31{
32 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N;
33};
34
38{
39public:
41 bool const is_axially_symmetric,
42 std::size_t n_local_size,
43 std::vector<unsigned>
44 dofIndex_to_localIndex)
45 : _element(element),
46 _is_axially_symmetric(is_axially_symmetric),
47 _dofIndex_to_localIndex(std::move(dofIndex_to_localIndex))
48 {
49 _local_u.resize(n_local_size);
50 _local_u_prev.resize(n_local_size);
51 _local_b.resize(_local_u.size());
52 _local_J.resize(_local_u.size(), _local_u.size());
53 }
54
55 void assemble(double const /*t*/, double const /*dt*/,
56 std::vector<double> const& /*local_x*/,
57 std::vector<double> const& /*local_x_prev*/,
58 std::vector<double>& /*local_M_data*/,
59 std::vector<double>& /*local_K_data*/,
60 std::vector<double>& /*local_rhs_data*/) override
61 {
63 "HydroMechanicsLocalAssembler: assembly without jacobian is not "
64 "implemented.");
65 }
66
67 void assembleWithJacobian(double const t, double const dt,
68 std::vector<double> const& local_x_,
69 std::vector<double> const& local_x_prev_,
70 std::vector<double>& /*local_M_data*/,
71 std::vector<double>& /*local_K_data*/,
72 std::vector<double>& local_b_data,
73 std::vector<double>& local_Jac_data) override
74 {
75 auto const local_dof_size = local_x_.size();
76
77 _local_u.setZero();
78 for (unsigned i = 0; i < local_dof_size; i++)
79 {
80 _local_u[_dofIndex_to_localIndex[i]] = local_x_[i];
81 }
82 _local_u_prev.setZero();
83 for (unsigned i = 0; i < local_dof_size; i++)
84 {
85 _local_u_prev[_dofIndex_to_localIndex[i]] = local_x_prev_[i];
86 }
87 _local_b.setZero();
88 _local_J.setZero();
89
91 _local_J);
92
93 local_b_data.resize(local_dof_size);
94 for (unsigned i = 0; i < local_dof_size; i++)
95 {
96 local_b_data[i] = _local_b[_dofIndex_to_localIndex[i]];
97 }
98
99 local_Jac_data.resize(local_dof_size * local_dof_size);
100 for (unsigned i = 0; i < local_dof_size; i++)
101 {
102 for (unsigned j = 0; j < local_dof_size; j++)
103 {
104 local_Jac_data[i * local_dof_size + j] = _local_J(
106 }
107 }
108 }
109
110 void postTimestepConcrete(Eigen::VectorXd const& local_x_,
111 Eigen::VectorXd const& /*local_x_prev*/,
112 const double t, double const dt,
113 int const /*process_id*/) override
114 {
115 _local_u.setZero();
116 for (Eigen::Index i = 0; i < local_x_.rows(); i++)
117 {
118 _local_u[_dofIndex_to_localIndex[i]] = local_x_[i];
119 }
120
122 }
123
124 virtual std::vector<double> const& getIntPtSigma(
125 const double t,
126 std::vector<GlobalVector*> const& x,
127 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
128 std::vector<double>& cache) const = 0;
129
130 virtual std::vector<double> const& getIntPtEpsilon(
131 const double t,
132 std::vector<GlobalVector*> const& x,
133 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
134 std::vector<double>& cache) const = 0;
135
136 virtual std::vector<double> const& getIntPtDarcyVelocity(
137 const double t,
138 std::vector<GlobalVector*> const& x,
139 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
140 std::vector<double>& cache) const = 0;
141
142 virtual std::vector<double> const& getIntPtFractureVelocity(
143 const double t,
144 std::vector<GlobalVector*> const& x,
145 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
146 std::vector<double>& cache) const = 0;
147
148 virtual std::vector<double> const& getIntPtFractureStress(
149 const double t,
150 std::vector<GlobalVector*> const& x,
151 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
152 std::vector<double>& cache) const = 0;
153
154 virtual std::vector<double> const& getIntPtFractureAperture(
155 const double t,
156 std::vector<GlobalVector*> const& x,
157 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
158 std::vector<double>& cache) const = 0;
159
160 virtual std::vector<double> const& getIntPtFracturePermeability(
161 const double t,
162 std::vector<GlobalVector*> const& x,
163 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
164 std::vector<double>& cache) const = 0;
165
166protected:
168 double const t, double const dt, Eigen::VectorXd const& local_u,
169 Eigen::VectorXd const& local_u_prev, Eigen::VectorXd& local_b,
170 Eigen::MatrixXd& local_J) = 0;
171
173 double const t, double const dt, Eigen::VectorXd const& local_u) = 0;
174
177
178private:
179 Eigen::VectorXd _local_u;
180 Eigen::VectorXd _local_u_prev;
181 Eigen::VectorXd _local_b;
182 Eigen::MatrixXd _local_J;
183 // a vector for mapping the index in the local DoF vector to the index in
184 // the complete local solution vector which also include nodes where DoF are
185 // not assigned.
186 std::vector<unsigned> const _dofIndex_to_localIndex;
187};
188
189} // namespace HydroMechanics
190} // namespace LIE
191} // 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, std::size_t n_local_size, std::vector< unsigned > dofIndex_to_localIndex)
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x_, std::vector< double > const &local_x_prev_, std::vector< double > &, std::vector< double > &, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
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
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