OGS
LiquidFlowLocalAssembler.h
Go to the documentation of this file.
1
13#pragma once
14
15#include <Eigen/Core>
16#include <vector>
17
18#include "LiquidFlowData.h"
32
33namespace ProcessLib
34{
35namespace LiquidFlow
36{
37template <typename GlobalDimNodalMatrixType>
39{
40 explicit IntegrationPointData(GlobalDimNodalMatrixType const& dNdx_,
41 double const& integration_weight_)
42 : dNdx(dNdx_), integration_weight(integration_weight_)
43 {
44 }
45
46 GlobalDimNodalMatrixType const dNdx;
47 double const integration_weight;
48
50};
51
52const unsigned NUM_NODAL_DOF = 1;
53
57{
58public:
59 virtual std::vector<double> const& getIntPtDarcyVelocity(
60 const double t,
61 std::vector<GlobalVector*> const& x,
62 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
63 std::vector<double>& cache) const = 0;
64};
65
66template <typename ShapeFunction, int GlobalDim>
68{
71
73 ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>;
74
82
83public:
85 MeshLib::Element const& element,
86 std::size_t const /*local_matrix_size*/,
87 NumLib::GenericIntegrationMethod const& integration_method,
88 bool const is_axially_symmetric,
89 LiquidFlowData const& process_data)
90 : _element(element),
91 _integration_method(integration_method),
92 _process_data(process_data)
93 {
94 unsigned const n_integration_points =
96 _ip_data.reserve(n_integration_points);
97
98 auto const& shape_matrices =
100 GlobalDim>(element, is_axially_symmetric,
102
105
106 double const aperture_size =
107 (_element.getDimension() == 3u)
108 ? 1.0
109 : _process_data.aperture_size(0.0, pos)[0];
110
111 for (unsigned ip = 0; ip < n_integration_points; ip++)
112 {
113 _ip_data.emplace_back(
114 shape_matrices[ip].dNdx,
116 shape_matrices[ip].integralMeasure *
117 shape_matrices[ip].detJ * aperture_size);
118 }
119 }
120
121 void assemble(double const t, double const dt,
122 std::vector<double> const& local_x,
123 std::vector<double> const& /*local_x_prev*/,
124 std::vector<double>& local_M_data,
125 std::vector<double>& local_K_data,
126 std::vector<double>& local_b_data) override;
127
130 Eigen::Vector3d getFlux(MathLib::Point3d const& p_local_coords,
131 double const t,
132 std::vector<double> const& local_x) const override;
133
134 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
135 const unsigned integration_point) const override
136 {
137 auto const& N =
139 .NsHigherOrder<typename ShapeFunction::MeshElement>();
140
141 // assumes N is stored contiguously in memory
142 return Eigen::Map<const Eigen::RowVectorXd>(
143 N[integration_point].data(), N[integration_point].size());
144 }
145
146 std::vector<double> const& getIntPtDarcyVelocity(
147 const double t,
148 std::vector<GlobalVector*> const& x,
149 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
150 std::vector<double>& velocity_cache) const override;
151
152private:
154
156 std::vector<IntegrationPointData<GlobalDimNodalMatrixType>> _ip_data;
157
163 {
165 Eigen::Map<NodalMatrixType>& local_K,
166 Eigen::Map<NodalVectorType>& local_b,
168 GlobalDimMatrixType const& permeability, double const mu,
169 double const rho_L, GlobalDimVectorType const& specific_body_force,
170 bool const has_gravity);
171
172 static Eigen::Matrix<double, GlobalDim, 1> calculateVelocity(
173 Eigen::Map<const NodalVectorType> const& local_p,
175 GlobalDimMatrixType const& permeability, double const mu,
176 double const rho_L, GlobalDimVectorType const& specific_body_force,
177 bool const has_gravity);
178 };
179
185 {
187 Eigen::Map<NodalMatrixType>& local_K,
188 Eigen::Map<NodalVectorType>& local_b,
190 GlobalDimMatrixType const& permeability, double const mu,
191 double const rho_L, GlobalDimVectorType const& specific_body_force,
192 bool const has_gravity);
193
194 static Eigen::Matrix<double, GlobalDim, 1> calculateVelocity(
195 Eigen::Map<const NodalVectorType> const& local_p,
197 GlobalDimMatrixType const& permeability, double const mu,
198 double const rho_L, GlobalDimVectorType const& specific_body_force,
199 bool const has_gravity);
200 };
201
202 template <typename LaplacianGravityVelocityCalculator>
203 void assembleMatrixAndVector(double const t, double const dt,
204 std::vector<double> const& local_x,
205 std::vector<double>& local_M_data,
206 std::vector<double>& local_K_data,
207 std::vector<double>& local_b_data);
208
209 template <typename LaplacianGravityVelocityCalculator,
210 typename VelocityCacheType>
212 const double t, const double dt, std::vector<double> const& local_x,
214 VelocityCacheType& darcy_velocity_at_ips) const;
215
216 template <typename VelocityCacheType>
217 void computeDarcyVelocity(bool const is_scalar_permeability, const double t,
218 const double dt,
219 std::vector<double> const& local_x,
221 VelocityCacheType& darcy_velocity_at_ips) const;
222
224};
225
226} // namespace LiquidFlow
227} // namespace ProcessLib
228
double getWeight() const
virtual constexpr unsigned getDimension() const =0
Get dimension of the mesh element.
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
auto const & NsHigherOrder() const
void setElementID(std::size_t element_id)
virtual std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< double > &cache) const =0
Eigen::Vector3d getFlux(MathLib::Point3d const &p_local_coords, double const t, std::vector< double > const &local_x) const override
std::vector< IntegrationPointData< GlobalDimNodalMatrixType > > _ip_data
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
typename LocalAssemblerTraits::LocalMatrix NodalMatrixType
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
LiquidFlowLocalAssembler(MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, LiquidFlowData const &process_data)
typename ShapeMatricesType::NodalRowVectorType NodalRowVectorType
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
void computeProjectedDarcyVelocity(const double t, const double dt, std::vector< double > const &local_x, ParameterLib::SpatialPosition const &pos, VelocityCacheType &darcy_velocity_at_ips) const
typename LocalAssemblerTraits::LocalVector NodalVectorType
void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override
NumLib::GenericIntegrationMethod const & _integration_method
void assembleMatrixAndVector(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
void computeDarcyVelocity(bool const is_scalar_permeability, const double t, const double dt, std::vector< double > const &local_x, ParameterLib::SpatialPosition const &pos, VelocityCacheType &darcy_velocity_at_ips) const
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
typename ShapeMatricesType::GlobalDimNodalMatrixType GlobalDimNodalMatrixType
std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &velocity_cache) const override
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
IntegrationPointData(GlobalDimNodalMatrixType const &dNdx_, double const &integration_weight_)
NumLib::ShapeMatrixCache shape_matrix_cache
caches for each mesh element type the shape matrix
ParameterLib::Parameter< double > const & aperture_size
static void calculateLaplacianAndGravityTerm(Eigen::Map< NodalMatrixType > &local_K, Eigen::Map< NodalVectorType > &local_b, IntegrationPointData< GlobalDimNodalMatrixType > const &ip_data, GlobalDimMatrixType const &permeability, double const mu, double const rho_L, GlobalDimVectorType const &specific_body_force, bool const has_gravity)
static Eigen::Matrix< double, GlobalDim, 1 > calculateVelocity(Eigen::Map< const NodalVectorType > const &local_p, IntegrationPointData< GlobalDimNodalMatrixType > const &ip_data, GlobalDimMatrixType const &permeability, double const mu, double const rho_L, GlobalDimVectorType const &specific_body_force, bool const has_gravity)
static void calculateLaplacianAndGravityTerm(Eigen::Map< NodalMatrixType > &local_K, Eigen::Map< NodalVectorType > &local_b, IntegrationPointData< GlobalDimNodalMatrixType > const &ip_data, GlobalDimMatrixType const &permeability, double const mu, double const rho_L, GlobalDimVectorType const &specific_body_force, bool const has_gravity)
static Eigen::Matrix< double, GlobalDim, 1 > calculateVelocity(Eigen::Map< const NodalVectorType > const &local_p, IntegrationPointData< GlobalDimNodalMatrixType > const &ip_data, GlobalDimMatrixType const &permeability, double const mu, double const rho_L, GlobalDimVectorType const &specific_body_force, bool const has_gravity)
Matrix< NNodes *NodalDOF, NNodes *NodalDOF > LocalMatrix