OGS
TwoPhaseFlowWithPrhoLocalAssembler.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <vector>
14
26
27namespace ProcessLib
28{
29namespace TwoPhaseFlowWithPrho
30{
31template <typename NodalMatrixType>
33{
36 : mat_property(material_property_),
37 sw(1.0),
38 rho_m(0.0),
39 dsw_dpg(0.0),
40 dsw_drho(0.0),
41 drhom_dpg(0.0),
42 drhom_drho(0.0)
43 {
44 }
46 double sw;
47 double rho_m;
48 double dsw_dpg;
49 double dsw_drho;
50 double drhom_dpg;
51 double drhom_drho;
53
55 NodalMatrixType massOperator;
56 NodalMatrixType diffusionOperator;
57
59};
60const unsigned NUM_NODAL_DOF = 2;
61
65{
66public:
67 virtual std::vector<double> const& getIntPtSaturation(
68 const double t,
69 std::vector<GlobalVector*> const& x,
70 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
71 std::vector<double>& cache) const = 0;
72
73 virtual std::vector<double> const& getIntPtNonWettingPressure(
74 const double t,
75 std::vector<GlobalVector*> const& x,
76 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
77 std::vector<double>& cache) const = 0;
78};
79
80template <typename ShapeFunction, int GlobalDim>
83{
86
88 ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>;
89
96
97public:
99 MeshLib::Element const& element,
100 std::size_t const /*local_matrix_size*/,
101 NumLib::GenericIntegrationMethod const& integration_method,
102 bool const is_axially_symmetric,
103 TwoPhaseFlowWithPrhoProcessData const& process_data)
104 : _element(element),
105 _integration_method(integration_method),
107 NumLib::initShapeMatrices<ShapeFunction, ShapeMatricesType,
108 GlobalDim>(
109 element, is_axially_symmetric, _integration_method)),
110 _process_data(process_data),
112 std::vector<double>(_integration_method.getNumberOfPoints())),
114 std::vector<double>(_integration_method.getNumberOfPoints()))
115 {
116 unsigned const n_integration_points =
118 _ip_data.reserve(n_integration_points);
119 for (unsigned ip = 0; ip < n_integration_points; ip++)
120 {
121 _ip_data.emplace_back(*_process_data._material);
122 auto const& sm = _shape_matrices[ip];
123 _ip_data[ip].integration_weight =
124 sm.integralMeasure * sm.detJ *
126 _ip_data[ip].massOperator.setZero(ShapeFunction::NPOINTS,
127 ShapeFunction::NPOINTS);
128 _ip_data[ip].diffusionOperator.setZero(ShapeFunction::NPOINTS,
129 ShapeFunction::NPOINTS);
130 _ip_data[ip].massOperator.noalias() =
131 sm.N.transpose() * sm.N * _ip_data[ip].integration_weight;
132 _ip_data[ip].diffusionOperator.noalias() =
133 sm.dNdx.transpose() * sm.dNdx * _ip_data[ip].integration_weight;
134 }
135 }
136
137 void assemble(double const t, double const /*dt*/,
138 std::vector<double> const& local_x,
139 std::vector<double> const& local_x_prev,
140 std::vector<double>& local_M_data,
141 std::vector<double>& local_K_data,
142 std::vector<double>& local_b_data) override;
143
144 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
145 const unsigned integration_point) const override
146 {
147 auto const& N = _shape_matrices[integration_point].N;
148
149 // assumes N is stored contiguously in memory
150 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
151 }
152
153 std::vector<double> const& getIntPtSaturation(
154 const double /*t*/,
155 std::vector<GlobalVector*> const& /*x*/,
156 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
157 std::vector<double>& /*cache*/) const override
158 {
159 assert(!_saturation.empty());
160 return _saturation;
161 }
162
163 std::vector<double> const& getIntPtNonWettingPressure(
164 const double /*t*/,
165 std::vector<GlobalVector*> const& /*x*/,
166 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
167 std::vector<double>& /*cache*/) const override
168 {
169 assert(!_pressure_nonwetting.empty());
171 }
172
173private:
175
177 std::vector<ShapeMatrices, Eigen::aligned_allocator<ShapeMatrices>>
179
181 std::vector<IntegrationPointData<NodalMatrixType>,
182 Eigen::aligned_allocator<IntegrationPointData<NodalMatrixType>>>
184
185 std::vector<double> _saturation;
186 std::vector<double> _pressure_nonwetting;
187 static const int nonwet_pressure_coeff_index = 0;
188 static const int cap_pressure_coeff_index = 1;
189
190 static const int nonwet_pressure_matrix_index = 0;
191 static const int cap_pressure_matrix_index = ShapeFunction::NPOINTS;
192
193 static const int nonwet_pressure_size = ShapeFunction::NPOINTS;
194 static const int cap_pressure_size = ShapeFunction::NPOINTS;
195};
196
197} // namespace TwoPhaseFlowWithPrho
198} // namespace ProcessLib
199
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 & getIntPtNonWettingPressure(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
void assemble(double const t, double const, 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
std::vector< IntegrationPointData< NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalMatrixType > > > _ip_data
TwoPhaseFlowWithPrhoLocalAssembler(MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, TwoPhaseFlowWithPrhoProcessData const &process_data)
std::vector< double > const & getIntPtNonWettingPressure(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
std::vector< double > _pressure_nonwetting
used for secondary variable output
std::vector< double > const & getIntPtSaturation(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
IntegrationPointData(TwoPhaseFlowWithPrhoMaterialProperties &material_property_)
std::unique_ptr< TwoPhaseFlowWithPrhoMaterialProperties > _material
Matrix< NNodes *NodalDOF, NNodes *NodalDOF > LocalMatrix