OGS
TwoPhaseFlowWithPPLocalAssembler.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 <vector>
7
18
19namespace ProcessLib
20{
21namespace TwoPhaseFlowWithPP
22{
23template <typename NodalRowVectorType, typename GlobalDimNodalMatrixType,
24 typename NodalMatrixType>
26{
27 explicit IntegrationPointData(NodalRowVectorType N_,
28 GlobalDimNodalMatrixType dNdx_,
29 double const& integration_weight_)
30 : N(std::move(N_)),
31 dNdx(std::move(dNdx_)),
32 integration_weight(integration_weight_)
33
34 {
35 }
36 NodalRowVectorType const N;
37 GlobalDimNodalMatrixType const dNdx;
38 const double integration_weight;
39
41};
42const unsigned NUM_NODAL_DOF = 2;
43
47{
48public:
49 virtual std::vector<double> const& getIntPtSaturation(
50 const double t,
51 std::vector<GlobalVector*> const& x,
52 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
53 std::vector<double>& cache) const = 0;
54
55 virtual std::vector<double> const& getIntPtWetPressure(
56 const double t,
57 std::vector<GlobalVector*> const& x,
58 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
59 std::vector<double>& cache) const = 0;
60};
61
62template <typename ShapeFunction, int GlobalDim>
65{
68
70 ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>;
72
81
82public:
84 MeshLib::Element const& element,
85 std::size_t const /*local_matrix_size*/,
86 NumLib::GenericIntegrationMethod const& integration_method,
87 bool const is_axially_symmetric,
88 TwoPhaseFlowWithPPProcessData const& process_data)
89 : _element(element),
90 _integration_method(integration_method),
91 _process_data(process_data),
93 std::vector<double>(_integration_method.getNumberOfPoints())),
95 std::vector<double>(_integration_method.getNumberOfPoints()))
96 {
97 unsigned const n_integration_points =
98 _integration_method.getNumberOfPoints();
99 _ip_data.reserve(n_integration_points);
100 auto const shape_matrices =
102 GlobalDim>(element, is_axially_symmetric,
104 for (unsigned ip = 0; ip < n_integration_points; ip++)
105 {
106 auto const& sm = shape_matrices[ip];
107 _ip_data.emplace_back(
108 sm.N, sm.dNdx,
109 sm.integralMeasure * sm.detJ *
110 _integration_method.getWeightedPoint(ip).getWeight());
111 }
112 }
113
114 void setInitialConditionsConcrete(Eigen::VectorXd const local_x,
115 double const t,
116 int const process_id) override;
117
118 void assemble(double const t, double const dt,
119 std::vector<double> const& local_x,
120 std::vector<double> const& local_x_prev,
121 std::vector<double>& local_M_data,
122 std::vector<double>& local_K_data,
123 std::vector<double>& local_b_data) override;
124
125 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
126 const unsigned integration_point) const override
127 {
128 auto const& N = _ip_data[integration_point].N;
129
130 // assumes N is stored contiguously in memory
131 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
132 }
133
134 std::vector<double> const& getIntPtSaturation(
135 const double /*t*/,
136 std::vector<GlobalVector*> const& /*x*/,
137 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
138 std::vector<double>& /*cache*/) const override
139 {
140 assert(!_saturation.empty());
141 return _saturation;
142 }
143
144 std::vector<double> const& getIntPtWetPressure(
145 const double /*t*/,
146 std::vector<GlobalVector*> const& /*x*/,
147 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
148 std::vector<double>& /*cache*/) const override
149 {
150 assert(!_pressure_wet.empty());
151 return _pressure_wet;
152 }
153
154private:
156
158
160 std::vector<
163 Eigen::aligned_allocator<IntegrationPointData<
166
167 // output vector for wetting phase saturation with
168 // respect to each integration point
169 std::vector<double> _saturation;
170 // output vector for wetting phase pressure with respect
171 // to each integration point
172 std::vector<double> _pressure_wet;
173 static const int nonwet_pressure_coeff_index = 0;
174 static const int cap_pressure_coeff_index = 1;
175
176 static const int nonwet_pressure_matrix_index = 0;
177 static const int cap_pressure_matrix_index = ShapeFunction::NPOINTS;
178
179 static const int nonwet_pressure_size = ShapeFunction::NPOINTS;
180 static const int cap_pressure_size = ShapeFunction::NPOINTS;
181};
182
183} // namespace TwoPhaseFlowWithPP
184} // namespace ProcessLib
185
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
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)
ProcessLib::LocalAssemblerTraits< ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim > LocalAssemblerTraits
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)
detail::LocalAssemblerTraitsFixed< ShpPol, NNodes, NodalDOF, Dim > LocalAssemblerTraits
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