OGS
RichardsComponentTransportFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <Eigen/Core>
14#include <vector>
15
28
29namespace ProcessLib
30{
31namespace RichardsComponentTransport
32{
33template <typename NodalRowVectorType, typename GlobalDimNodalMatrixType,
34 typename NodalMatrixType>
36{
37 IntegrationPointData(NodalRowVectorType const& N_,
38 GlobalDimNodalMatrixType const& dNdx_,
39 double const& integration_weight_,
40 NodalMatrixType const mass_operator_)
41 : N(N_),
42 dNdx(dNdx_),
43 integration_weight(integration_weight_),
44 mass_operator(mass_operator_)
45 {
46 }
47
48 NodalRowVectorType const N;
49 GlobalDimNodalMatrixType const dNdx;
50 double const integration_weight;
51 NodalMatrixType const mass_operator;
52
54};
55
56const unsigned NUM_NODAL_DOF = 2;
57
61{
62public:
63 virtual std::vector<double> const& getIntPtSaturation(
64 const double t,
65 std::vector<GlobalVector*> const& x,
66 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
67 std::vector<double>& cache) const = 0;
68
69 virtual std::vector<double> const& getIntPtDarcyVelocity(
70 const double t,
71 std::vector<GlobalVector*> const& x,
72 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
73 std::vector<double>& cache) const = 0;
74};
75
76template <typename ShapeFunction, int GlobalDim>
79{
82
83 using LocalMatrixType = typename ShapeMatricesType::template MatrixType<
84 NUM_NODAL_DOF * ShapeFunction::NPOINTS,
85 NUM_NODAL_DOF * ShapeFunction::NPOINTS>;
87 typename ShapeMatricesType::template VectorType<NUM_NODAL_DOF *
88 ShapeFunction::NPOINTS>;
89
92
98
99public:
101 MeshLib::Element const& element,
102 std::size_t const local_matrix_size,
103 NumLib::GenericIntegrationMethod const& integration_method,
104 bool is_axially_symmetric,
105 RichardsComponentTransportProcessData const& process_data,
106 ProcessVariable const& transport_process_variable);
107
108 void assemble(double const t, double const dt,
109 std::vector<double> const& local_x,
110 std::vector<double> const& local_x_prev,
111 std::vector<double>& local_M_data,
112 std::vector<double>& local_K_data,
113 std::vector<double>& local_b_data) override;
114
115 std::vector<double> const& getIntPtDarcyVelocity(
116 const double t,
117 std::vector<GlobalVector*> const& x,
118 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
119 std::vector<double>& cache) const override;
120
121 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
122 const unsigned integration_point) const override;
123
124 std::vector<double> const& getIntPtSaturation(
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 override;
129
130private:
131 unsigned const _element_id;
133
136
137 std::vector<
140 Eigen::aligned_allocator<IntegrationPointData<
143
144 static const int concentration_index = 0;
145 static const int concentration_size = ShapeFunction::NPOINTS;
146 static const int pressure_index = ShapeFunction::NPOINTS;
147 static const int pressure_size = ShapeFunction::NPOINTS;
148};
149
150} // namespace RichardsComponentTransport
151} // namespace ProcessLib
152
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
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_)