OGS
TwoPhaseFlowWithPPLocalAssembler.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <vector>
14
25
26namespace ProcessLib
27{
28namespace TwoPhaseFlowWithPP
29{
30template <typename NodalRowVectorType, typename GlobalDimNodalMatrixType,
31 typename NodalMatrixType>
33{
34 explicit IntegrationPointData(NodalRowVectorType N_,
35 GlobalDimNodalMatrixType dNdx_,
36 double const& integration_weight_,
37 NodalMatrixType const massOperator_)
38 : N(std::move(N_)),
39 dNdx(std::move(dNdx_)),
40 integration_weight(integration_weight_),
41 massOperator(massOperator_)
42
43 {
44 }
45 NodalRowVectorType const N;
46 GlobalDimNodalMatrixType const dNdx;
47 const double integration_weight;
48 NodalMatrixType const massOperator;
49
51};
52const unsigned NUM_NODAL_DOF = 2;
53
57{
58public:
59 virtual std::vector<double> const& getIntPtSaturation(
60 const double t,
61 std::vector<GlobalVector*> const& x,
62 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
63 std::vector<double>& cache) const = 0;
64
65 virtual std::vector<double> const& getIntPtWetPressure(
66 const double t,
67 std::vector<GlobalVector*> const& x,
68 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
69 std::vector<double>& cache) const = 0;
70};
71
72template <typename ShapeFunction, int GlobalDim>
75{
78
80 ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>;
82
91
92public:
94 MeshLib::Element const& element,
95 std::size_t const /*local_matrix_size*/,
96 NumLib::GenericIntegrationMethod const& integration_method,
97 bool const is_axially_symmetric,
98 TwoPhaseFlowWithPPProcessData const& process_data)
99 : _element(element),
100 _integration_method(integration_method),
101 _process_data(process_data),
103 std::vector<double>(_integration_method.getNumberOfPoints())),
105 std::vector<double>(_integration_method.getNumberOfPoints()))
106 {
107 unsigned const n_integration_points =
109 _ip_data.reserve(n_integration_points);
110 auto const shape_matrices =
112 GlobalDim>(element, is_axially_symmetric,
114 for (unsigned ip = 0; ip < n_integration_points; ip++)
115 {
116 auto const& sm = shape_matrices[ip];
117 _ip_data.emplace_back(
118 sm.N, sm.dNdx,
119 sm.integralMeasure * sm.detJ *
121 sm.N.transpose() * sm.N * sm.integralMeasure * sm.detJ *
123 }
124 }
125
126 void assemble(double const t, double const dt,
127 std::vector<double> const& local_x,
128 std::vector<double> const& local_x_prev,
129 std::vector<double>& local_M_data,
130 std::vector<double>& local_K_data,
131 std::vector<double>& local_b_data) override;
132
133 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
134 const unsigned integration_point) const override
135 {
136 auto const& N = _ip_data[integration_point].N;
137
138 // assumes N is stored contiguously in memory
139 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
140 }
141
142 std::vector<double> const& getIntPtSaturation(
143 const double /*t*/,
144 std::vector<GlobalVector*> const& /*x*/,
145 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
146 std::vector<double>& /*cache*/) const override
147 {
148 assert(!_saturation.empty());
149 return _saturation;
150 }
151
152 std::vector<double> const& getIntPtWetPressure(
153 const double /*t*/,
154 std::vector<GlobalVector*> const& /*x*/,
155 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
156 std::vector<double>& /*cache*/) const override
157 {
158 assert(!_pressure_wet.empty());
159 return _pressure_wet;
160 }
161
162private:
164
166
168 std::vector<
171 Eigen::aligned_allocator<IntegrationPointData<
174
175 // output vector for wetting phase saturation with
176 // respect to each integration point
177 std::vector<double> _saturation;
178 // output vector for wetting phase pressure with respect
179 // to each integration point
180 std::vector<double> _pressure_wet;
181 static const int nonwet_pressure_coeff_index = 0;
182 static const int cap_pressure_coeff_index = 1;
183
184 static const int nonwet_pressure_matrix_index = 0;
185 static const int cap_pressure_matrix_index = ShapeFunction::NPOINTS;
186
187 static const int nonwet_pressure_size = ShapeFunction::NPOINTS;
188 static const int cap_pressure_size = ShapeFunction::NPOINTS;
189};
190
191} // namespace TwoPhaseFlowWithPP
192} // namespace ProcessLib
193
double getWeight() const
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
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
virtual std::vector< double > const & getIntPtWetPressure(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
std::vector< double > const & getIntPtSaturation(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
TwoPhaseFlowWithPPLocalAssembler(MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, TwoPhaseFlowWithPPProcessData const &process_data)
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType > > > _ip_data
typename ShapeMatricesType::GlobalDimNodalMatrixType GlobalDimNodalMatrixType
std::vector< double > const & getIntPtWetPressure(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) 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
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
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< 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 N_, GlobalDimNodalMatrixType dNdx_, double const &integration_weight_, NodalMatrixType const massOperator_)
Matrix< NNodes *NodalDOF, NNodes *NodalDOF > LocalMatrix