OGS
LiquidFlowLocalAssembler.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 <Eigen/Core>
7#include <vector>
8
9#include "LiquidFlowData.h"
23
24namespace ProcessLib
25{
26namespace LiquidFlow
27{
28template <typename GlobalDimNodalMatrixType>
30{
31 explicit IntegrationPointData(GlobalDimNodalMatrixType const& dNdx_,
32 double const& integration_weight_)
33 : dNdx(dNdx_), integration_weight(integration_weight_)
34 {
35 }
36
37 GlobalDimNodalMatrixType const dNdx;
38 double const integration_weight;
39
41};
42
43const unsigned NUM_NODAL_DOF = 1;
44
48{
49public:
50 virtual std::vector<double> const& getIntPtDarcyVelocity(
51 const double t,
52 std::vector<GlobalVector*> const& x,
53 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
54 std::vector<double>& cache) const = 0;
55};
56
57template <typename ShapeFunction, int GlobalDim>
59{
62
64 ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>;
65
73
74public:
76 MeshLib::Element const& element,
77 std::size_t const /*local_matrix_size*/,
78 NumLib::GenericIntegrationMethod const& integration_method,
79 bool const is_axially_symmetric,
80 LiquidFlowData const& process_data)
81 : _element(element),
82 _integration_method(integration_method),
83 _process_data(process_data)
84 {
85 unsigned const n_integration_points =
86 _integration_method.getNumberOfPoints();
87 _ip_data.reserve(n_integration_points);
88
89 auto const& shape_matrices =
91 GlobalDim>(element, is_axially_symmetric,
93
95 pos.setElementID(_element.getID());
96
97 double const aperture_size =
98 (_element.getDimension() == 3u)
99 ? 1.0
100 : _process_data.aperture_size(0.0, pos)[0];
101
102 for (unsigned ip = 0; ip < n_integration_points; ip++)
103 {
104 _ip_data.emplace_back(
105 shape_matrices[ip].dNdx,
106 _integration_method.getWeightedPoint(ip).getWeight() *
107 shape_matrices[ip].integralMeasure *
108 shape_matrices[ip].detJ * aperture_size);
109 }
110 }
111
112 void assemble(double const t, double const dt,
113 std::vector<double> const& local_x,
114 std::vector<double> const& /*local_x_prev*/,
115 std::vector<double>& local_M_data,
116 std::vector<double>& local_K_data,
117 std::vector<double>& local_b_data) override;
118
121 Eigen::Vector3d getFlux(MathLib::Point3d const& p_local_coords,
122 double const t,
123 std::vector<double> const& local_x) const override;
124
125 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
126 const unsigned integration_point) const override
127 {
128 auto const& N =
129 _process_data.shape_matrix_cache
130 .NsHigherOrder<typename ShapeFunction::MeshElement>();
131
132 // assumes N is stored contiguously in memory
133 return Eigen::Map<const Eigen::RowVectorXd>(
134 N[integration_point].data(), N[integration_point].size());
135 }
136
137 std::vector<double> const& getIntPtDarcyVelocity(
138 const double t,
139 std::vector<GlobalVector*> const& x,
140 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
141 std::vector<double>& velocity_cache) const override;
142
143private:
145
147 std::vector<IntegrationPointData<GlobalDimNodalMatrixType>> _ip_data;
148
154 {
156 Eigen::Map<NodalMatrixType>& local_K,
157 Eigen::Map<NodalVectorType>& local_b,
159 GlobalDimMatrixType const& permeability_with_density_factor,
160 double const mu, double const rho_L,
161 GlobalDimVectorType const& specific_body_force,
162 bool const has_gravity);
163
164 static Eigen::Matrix<double, GlobalDim, 1> calculateVelocity(
165 Eigen::Map<const NodalVectorType> const& local_p,
167 GlobalDimMatrixType const& permeability, double const mu,
168 double const rho_L, GlobalDimVectorType const& specific_body_force,
169 bool const has_gravity);
170 };
171
177 {
179 Eigen::Map<NodalMatrixType>& local_K,
180 Eigen::Map<NodalVectorType>& local_b,
182 GlobalDimMatrixType const& permeability_with_density_factor,
183 double const mu, double const rho_L,
184 GlobalDimVectorType const& specific_body_force,
185 bool const has_gravity);
186
187 static Eigen::Matrix<double, GlobalDim, 1> calculateVelocity(
188 Eigen::Map<const NodalVectorType> const& local_p,
190 GlobalDimMatrixType const& permeability, double const mu,
191 double const rho_L, GlobalDimVectorType const& specific_body_force,
192 bool const has_gravity);
193 };
194
195 template <typename LaplacianGravityVelocityCalculator>
196 void assembleMatrixAndVector(double const t, double const dt,
197 std::vector<double> const& local_x,
198 std::vector<double>& local_M_data,
199 std::vector<double>& local_K_data,
200 std::vector<double>& local_b_data);
201
202 template <typename LaplacianGravityVelocityCalculator,
203 typename VelocityCacheType>
205 const double t, const double dt, std::vector<double> const& local_x,
207 VelocityCacheType& darcy_velocity_at_ips) const;
208
209 template <typename VelocityCacheType>
210 void computeDarcyVelocity(bool const is_scalar_permeability, const double t,
211 const double dt,
212 std::vector<double> const& local_x,
214 VelocityCacheType& darcy_velocity_at_ips) const;
215
217};
218
219} // namespace LiquidFlow
220} // namespace ProcessLib
221
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
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
ProcessLib::LocalAssemblerTraits< ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim > LocalAssemblerTraits
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)
detail::LocalAssemblerTraitsFixed< ShpPol, NNodes, NodalDOF, Dim > LocalAssemblerTraits
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_)
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