OGS
ThermalTwoPhaseFlowWithPPLocalAssembler.h
Go to the documentation of this file.
1 
11 #pragma once
12 
13 #include <vector>
14 
21 #include "ParameterLib/Parameter.h"
25 
26 namespace ProcessLib
27 {
28 namespace ThermalTwoPhaseFlowWithPP
29 {
30 template <typename NodalRowVectorType, typename GlobalDimNodalMatrixType,
31  typename NodalMatrixType>
33 {
35  NodalRowVectorType const& N_, GlobalDimNodalMatrixType const& dNdx_,
36  ThermalTwoPhaseFlowWithPPMaterialProperties const& material_property_,
37  double const& integration_weight_, NodalMatrixType const mass_operator_,
38  NodalMatrixType const diffusion_operator_)
39  : N(N_),
40  dNdx(dNdx_),
41  mat_property(material_property_),
42  integration_weight(integration_weight_),
43  mass_operator(mass_operator_),
44  diffusion_operator(diffusion_operator_)
45  {
46  }
47  NodalRowVectorType const N;
48  GlobalDimNodalMatrixType const dNdx;
50  double const integration_weight;
51  NodalMatrixType const mass_operator;
52  NodalMatrixType const diffusion_operator;
53 
55 };
56 const unsigned NUM_NODAL_DOF = 3;
57 
61 {
62 public:
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& getIntPtWettingPressure(
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 
76 template <typename ShapeFunction, typename IntegrationMethod, int GlobalDim>
79 {
82 
84  ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>;
86 
95 
96 public:
98  MeshLib::Element const& element,
99  std::size_t const /*local_matrix_size*/,
100  bool const is_axially_symmetric,
101  unsigned const integration_order,
102  ThermalTwoPhaseFlowWithPPProcessData const& process_data)
103  : _element(element),
104  _integration_method(integration_order),
105  _process_data(process_data),
106  _saturation(
107  std::vector<double>(_integration_method.getNumberOfPoints())),
109  std::vector<double>(_integration_method.getNumberOfPoints()))
110  {
111  unsigned const n_integration_points =
112  _integration_method.getNumberOfPoints();
113  _ip_data.reserve(n_integration_points);
114  auto const shape_matrices =
116  GlobalDim>(element, is_axially_symmetric,
118  for (unsigned ip = 0; ip < n_integration_points; ip++)
119  {
120  auto const& sm = shape_matrices[ip];
121  const double integration_factor = sm.integralMeasure * sm.detJ;
122  _ip_data.emplace_back(
123  sm.N, sm.dNdx, *_process_data.material,
124  sm.integralMeasure * sm.detJ *
125  _integration_method.getWeightedPoint(ip).getWeight(),
126  sm.N.transpose() * sm.N * integration_factor *
127  _integration_method.getWeightedPoint(ip).getWeight(),
128  sm.dNdx.transpose() * sm.dNdx * integration_factor *
129  _integration_method.getWeightedPoint(ip).getWeight());
130  }
131  }
132 
133  void assemble(double const t, double const dt,
134  std::vector<double> const& local_x,
135  std::vector<double> const& local_xdot,
136  std::vector<double>& local_M_data,
137  std::vector<double>& local_K_data,
138  std::vector<double>& local_b_data) override;
139 
140  Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
141  const unsigned integration_point) const override
142  {
143  auto const& N = _ip_data[integration_point].N;
144 
145  // assumes N is stored contiguously in memory
146  return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
147  }
148 
149  std::vector<double> const& getIntPtSaturation(
150  const double /*t*/,
151  std::vector<GlobalVector*> const& /*x*/,
152  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
153  std::vector<double>& /*cache*/) const override
154  {
155  assert(!_saturation.empty());
156  return _saturation;
157  }
158 
159  std::vector<double> const& getIntPtWettingPressure(
160  const double /*t*/,
161  std::vector<GlobalVector*> const& /*x*/,
162  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
163  std::vector<double>& /*cache*/) const override
164  {
165  assert(!_pressure_wetting.empty());
166  return _pressure_wetting;
167  }
168 
169 private:
171 
172  IntegrationMethod const _integration_method;
173 
175  std::vector<
178  Eigen::aligned_allocator<IntegrationPointData<
181 
182  std::vector<double> _saturation;
183  std::vector<double> _pressure_wetting;
184 
185  static const int nonwet_pressure_matrix_index = 0;
186  static const int cap_pressure_matrix_index = ShapeFunction::NPOINTS;
187  static const int temperature_matrix_index = 2 * ShapeFunction::NPOINTS;
188 
189  static const int nonwet_pressure_size = ShapeFunction::NPOINTS;
190  static const int cap_pressure_size = ShapeFunction::NPOINTS;
191  static const int temperature_size = ShapeFunction::NPOINTS;
192 };
193 
194 } // namespace ThermalTwoPhaseFlowWithPP
195 } // namespace ProcessLib
196 
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 & getIntPtWettingPressure(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
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType > > > _ip_data
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 dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override
std::vector< double > const & getIntPtWettingPressure(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
ThermalTwoPhaseFlowWithPPLocalAssembler(MeshLib::Element const &element, std::size_t const, bool const is_axially_symmetric, unsigned const integration_order, ThermalTwoPhaseFlowWithPPProcessData const &process_data)
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)
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
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_, ThermalTwoPhaseFlowWithPPMaterialProperties const &material_property_, double const &integration_weight_, NodalMatrixType const mass_operator_, NodalMatrixType const diffusion_operator_)
Matrix< NNodes *NodalDOF, NNodes *NodalDOF > LocalMatrix