OGS
RichardsComponentTransportFEM.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
21
22namespace ProcessLib
23{
25{
26template <typename NodalRowVectorType, typename GlobalDimNodalMatrixType,
27 typename NodalMatrixType>
29{
30 IntegrationPointData(NodalRowVectorType const& N_,
31 GlobalDimNodalMatrixType const& dNdx_,
32 double const& integration_weight_,
33 NodalMatrixType const mass_operator_)
34 : N(N_),
35 dNdx(dNdx_),
36 integration_weight(integration_weight_),
37 mass_operator(mass_operator_)
38 {
39 }
40
41 NodalRowVectorType const N;
42 GlobalDimNodalMatrixType const dNdx;
43 double const integration_weight;
44 NodalMatrixType const mass_operator;
45
47};
48
49const unsigned NUM_NODAL_DOF = 2;
50
54{
55public:
56 virtual std::vector<double> const& getIntPtSaturation(
57 const double t,
58 std::vector<GlobalVector*> const& x,
59 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
60 std::vector<double>& cache) const = 0;
61
62 virtual std::vector<double> const& getIntPtDarcyVelocity(
63 const double t,
64 std::vector<GlobalVector*> const& x,
65 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
66 std::vector<double>& cache) const = 0;
67};
68
69template <typename ShapeFunction, int GlobalDim>
72{
75
76 using LocalMatrixType = typename ShapeMatricesType::template MatrixType<
77 NUM_NODAL_DOF * ShapeFunction::NPOINTS,
78 NUM_NODAL_DOF * ShapeFunction::NPOINTS>;
80 typename ShapeMatricesType::template VectorType<NUM_NODAL_DOF *
81 ShapeFunction::NPOINTS>;
82
85
91
92public:
94 MeshLib::Element const& element,
95 std::size_t const local_matrix_size,
96 NumLib::GenericIntegrationMethod const& integration_method,
97 bool is_axially_symmetric,
98 RichardsComponentTransportProcessData const& process_data,
99 ProcessVariable const& transport_process_variable);
100
101 void assemble(double const t, double const dt,
102 std::vector<double> const& local_x,
103 std::vector<double> const& local_x_prev,
104 std::vector<double>& local_M_data,
105 std::vector<double>& local_K_data,
106 std::vector<double>& local_b_data) override;
107
108 std::vector<double> const& getIntPtDarcyVelocity(
109 const double t,
110 std::vector<GlobalVector*> const& x,
111 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
112 std::vector<double>& cache) const override;
113
114 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
115 const unsigned integration_point) const override;
116
117 std::vector<double> const& getIntPtSaturation(
118 const double t,
119 std::vector<GlobalVector*> const& x,
120 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
121 std::vector<double>& cache) const override;
122
123private:
126
129
130 std::vector<
133 Eigen::aligned_allocator<IntegrationPointData<
136
137 static const int concentration_index = 0;
138 static const int concentration_size = ShapeFunction::NPOINTS;
139 static const int pressure_index = ShapeFunction::NPOINTS;
140 static const int pressure_size = ShapeFunction::NPOINTS;
141};
142
143} // namespace RichardsComponentTransport
144} // namespace ProcessLib
145
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType > > > _ip_data
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
typename ShapeMatricesType::GlobalDimNodalMatrixType GlobalDimNodalMatrixType
typename ShapeMatricesType::template VectorType< NUM_NODAL_DOF * ShapeFunction::NPOINTS > LocalVectorType
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
LocalAssemblerData(MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool is_axially_symmetric, RichardsComponentTransportProcessData const &process_data, ProcessVariable const &transport_process_variable)
std::vector< double > const & getIntPtSaturation(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
typename ShapeMatricesType::NodalRowVectorType NodalRowVectorType
typename ShapeMatricesType::template MatrixType< NUM_NODAL_DOF *ShapeFunction::NPOINTS, NUM_NODAL_DOF *ShapeFunction::NPOINTS > LocalMatrixType
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 override
void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
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 & getIntPtSaturation(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
IntegrationPointData(NodalRowVectorType const &N_, GlobalDimNodalMatrixType const &dNdx_, double const &integration_weight_, NodalMatrixType const mass_operator_)