OGS
LiquidFlowLocalAssembler.h
Go to the documentation of this file.
1
13#pragma once
14
15#include <Eigen/Core>
16#include <tuple>
17#include <vector>
18
19#include "LiquidFlowData.h"
33
34namespace ProcessLib
35{
36namespace LiquidFlow
37{
38template <typename GlobalDimNodalMatrixType>
40{
41 explicit IntegrationPointData(GlobalDimNodalMatrixType const& dNdx_,
42 double const& integration_weight_)
43 : dNdx(dNdx_), integration_weight(integration_weight_)
44 {
45 }
46
47 GlobalDimNodalMatrixType const dNdx;
48 double const integration_weight;
49
51};
52
53const unsigned NUM_NODAL_DOF = 1;
54
58{
59public:
60 virtual std::vector<double> const& getIntPtDarcyVelocity(
61 const double t,
62 std::vector<GlobalVector*> const& x,
63 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
64 std::vector<double>& cache) const = 0;
65};
66
67template <typename ShapeFunction, int GlobalDim>
69{
72
74 ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>;
75
83
84public:
86 MeshLib::Element const& element,
87 std::size_t const /*local_matrix_size*/,
88 NumLib::GenericIntegrationMethod const& integration_method,
89 bool const is_axially_symmetric,
90 LiquidFlowData const& process_data)
91 : _element(element),
92 _integration_method(integration_method),
93 _process_data(process_data)
94 {
95 unsigned const n_integration_points =
97 _ip_data.reserve(n_integration_points);
98
99 auto const& shape_matrices =
101 GlobalDim>(element, is_axially_symmetric,
103
106
107 double const aperture_size =
108 (_element.getDimension() == 3u)
109 ? 1.0
110 : _process_data.aperture_size(0.0, pos)[0];
111
112 for (unsigned ip = 0; ip < n_integration_points; ip++)
113 {
114 _ip_data.emplace_back(
115 shape_matrices[ip].dNdx,
117 shape_matrices[ip].integralMeasure *
118 shape_matrices[ip].detJ * aperture_size);
119 }
120 }
121
122 void assemble(double const t, double const dt,
123 std::vector<double> const& local_x,
124 std::vector<double> const& /*local_x_prev*/,
125 std::vector<double>& local_M_data,
126 std::vector<double>& local_K_data,
127 std::vector<double>& local_b_data) override;
128
131 Eigen::Vector3d getFlux(MathLib::Point3d const& p_local_coords,
132 double const t,
133 std::vector<double> const& local_x) const override;
134
135 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
136 const unsigned integration_point) const override
137 {
138 auto const& N =
140 .NsHigherOrder<typename ShapeFunction::MeshElement>();
141
142 // assumes N is stored contiguously in memory
143 return Eigen::Map<const Eigen::RowVectorXd>(
144 N[integration_point].data(), N[integration_point].size());
145 }
146
147 std::vector<double> const& getIntPtDarcyVelocity(
148 const double t,
149 std::vector<GlobalVector*> const& x,
150 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
151 std::vector<double>& velocity_cache) const override;
152
153private:
155
157 std::vector<IntegrationPointData<GlobalDimNodalMatrixType>> _ip_data;
158
164 {
166 Eigen::Map<NodalMatrixType>& local_K,
167 Eigen::Map<NodalVectorType>& local_b,
169 GlobalDimMatrixType const& permeability_with_density_factor,
170 double const mu, double const rho_L,
171 GlobalDimVectorType const& specific_body_force,
172 bool const has_gravity);
173
174 static Eigen::Matrix<double, GlobalDim, 1> calculateVelocity(
175 Eigen::Map<const NodalVectorType> const& local_p,
177 GlobalDimMatrixType const& permeability, double const mu,
178 double const rho_L, GlobalDimVectorType const& specific_body_force,
179 bool const has_gravity);
180 };
181
187 {
189 Eigen::Map<NodalMatrixType>& local_K,
190 Eigen::Map<NodalVectorType>& local_b,
192 GlobalDimMatrixType const& permeability_with_density_factor,
193 double const mu, double const rho_L,
194 GlobalDimVectorType const& specific_body_force,
195 bool const has_gravity);
196
197 static Eigen::Matrix<double, GlobalDim, 1> calculateVelocity(
198 Eigen::Map<const NodalVectorType> const& local_p,
200 GlobalDimMatrixType const& permeability, double const mu,
201 double const rho_L, GlobalDimVectorType const& specific_body_force,
202 bool const has_gravity);
203 };
204
205 template <typename LaplacianGravityVelocityCalculator>
206 void assembleMatrixAndVector(double const t, double const dt,
207 std::vector<double> const& local_x,
208 std::vector<double>& local_M_data,
209 std::vector<double>& local_K_data,
210 std::vector<double>& local_b_data);
211
212 template <typename LaplacianGravityVelocityCalculator,
213 typename VelocityCacheType>
215 const double t, const double dt, std::vector<double> const& local_x,
217 VelocityCacheType& darcy_velocity_at_ips) const;
218
219 template <typename VelocityCacheType>
220 void computeDarcyVelocity(bool const is_scalar_permeability, const double t,
221 const double dt,
222 std::vector<double> const& local_x,
224 VelocityCacheType& darcy_velocity_at_ips) const;
225
227};
228
229std::tuple<double, double> getFluidDensityAndViscosity(
230 double const t, double const dt, ParameterLib::SpatialPosition const& pos,
231 MaterialPropertyLib::Phase const& fluid_phase,
233
234} // namespace LiquidFlow
235} // namespace ProcessLib
236
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)
std::tuple< double, double > getFluidDensityAndViscosity(double const t, double const dt, ParameterLib::SpatialPosition const &pos, MaterialPropertyLib::Phase const &fluid_phase, MaterialPropertyLib::VariableArray &vars)
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_with_density_factor, 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_with_density_factor, 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