Loading [MathJax]/extensions/tex2jax.js
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 : N(std::move(N_)),
38 dNdx(std::move(dNdx_)),
39 integration_weight(integration_weight_)
40
41 {
42 }
43 NodalRowVectorType const N;
44 GlobalDimNodalMatrixType const dNdx;
45 const double integration_weight;
46
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& getIntPtWetPressure(
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
77 ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>;
79
88
89public:
91 MeshLib::Element const& element,
92 std::size_t const /*local_matrix_size*/,
93 NumLib::GenericIntegrationMethod const& integration_method,
94 bool const is_axially_symmetric,
95 TwoPhaseFlowWithPPProcessData const& process_data)
96 : _element(element),
97 _integration_method(integration_method),
98 _process_data(process_data),
100 std::vector<double>(_integration_method.getNumberOfPoints())),
102 std::vector<double>(_integration_method.getNumberOfPoints()))
103 {
104 unsigned const n_integration_points =
106 _ip_data.reserve(n_integration_points);
107 auto const shape_matrices =
109 GlobalDim>(element, is_axially_symmetric,
111 for (unsigned ip = 0; ip < n_integration_points; ip++)
112 {
113 auto const& sm = shape_matrices[ip];
114 _ip_data.emplace_back(
115 sm.N, sm.dNdx,
116 sm.integralMeasure * sm.detJ *
118 }
119 }
120
121 void setInitialConditionsConcrete(Eigen::VectorXd const local_x,
122 double const t,
123 int const process_id) override;
124
125 void assemble(double const t, double const dt,
126 std::vector<double> const& local_x,
127 std::vector<double> const& local_x_prev,
128 std::vector<double>& local_M_data,
129 std::vector<double>& local_K_data,
130 std::vector<double>& local_b_data) override;
131
132 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
133 const unsigned integration_point) const override
134 {
135 auto const& N = _ip_data[integration_point].N;
136
137 // assumes N is stored contiguously in memory
138 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
139 }
140
141 std::vector<double> const& getIntPtSaturation(
142 const double /*t*/,
143 std::vector<GlobalVector*> const& /*x*/,
144 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
145 std::vector<double>& /*cache*/) const override
146 {
147 assert(!_saturation.empty());
148 return _saturation;
149 }
150
151 std::vector<double> const& getIntPtWetPressure(
152 const double /*t*/,
153 std::vector<GlobalVector*> const& /*x*/,
154 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
155 std::vector<double>& /*cache*/) const override
156 {
157 assert(!_pressure_wet.empty());
158 return _pressure_wet;
159 }
160
161private:
163
165
167 std::vector<
170 Eigen::aligned_allocator<IntegrationPointData<
173
174 // output vector for wetting phase saturation with
175 // respect to each integration point
176 std::vector<double> _saturation;
177 // output vector for wetting phase pressure with respect
178 // to each integration point
179 std::vector<double> _pressure_wet;
180 static const int nonwet_pressure_coeff_index = 0;
181 static const int cap_pressure_coeff_index = 1;
182
183 static const int nonwet_pressure_matrix_index = 0;
184 static const int cap_pressure_matrix_index = ShapeFunction::NPOINTS;
185
186 static const int nonwet_pressure_size = ShapeFunction::NPOINTS;
187 static const int cap_pressure_size = ShapeFunction::NPOINTS;
188};
189
190} // namespace TwoPhaseFlowWithPP
191} // namespace ProcessLib
192
constexpr 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
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
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_)
Matrix< NNodes *NodalDOF, NNodes *NodalDOF > LocalMatrix