29namespace ThermalTwoPhaseFlowWithPP
31template <
typename NodalRowVectorType,
typename GlobalDimNodalMatrixType,
32 typename NodalMatrixType>
36 GlobalDimNodalMatrixType
const& dNdx_,
37 double const& integration_weight_,
38 NodalMatrixType
const mass_operator_,
39 NodalMatrixType
const diffusion_operator_)
47 NodalRowVectorType
const N;
48 GlobalDimNodalMatrixType
const dNdx;
64 std::vector<GlobalVector*>
const& x,
65 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
66 std::vector<double>& cache)
const = 0;
70 std::vector<GlobalVector*>
const& x,
71 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
72 std::vector<double>& cache)
const = 0;
76 std::vector<GlobalVector*>
const& x,
77 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
78 std::vector<double>& cache)
const = 0;
82 std::vector<GlobalVector*>
const& x,
83 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
84 std::vector<double>& cache)
const = 0;
88 std::vector<GlobalVector*>
const& x,
89 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
90 std::vector<double>& cache)
const = 0;
93template <
typename ShapeFunction,
int GlobalDim>
118 bool const is_axially_symmetric,
134 unsigned const n_integration_points =
136 _ip_data.reserve(n_integration_points);
137 auto const shape_matrices =
139 GlobalDim>(element, is_axially_symmetric,
141 for (
unsigned ip = 0; ip < n_integration_points; ip++)
143 auto const& sm = shape_matrices[ip];
144 const double integration_factor = sm.integralMeasure * sm.detJ;
147 sm.integralMeasure * sm.detJ *
149 sm.N.transpose() * sm.N * integration_factor *
151 sm.dNdx.transpose() * sm.dNdx * integration_factor *
156 void assemble(
double const t,
double const dt,
157 std::vector<double>
const& local_x,
158 std::vector<double>
const& local_x_prev,
159 std::vector<double>& local_M_data,
160 std::vector<double>& local_K_data,
161 std::vector<double>& local_b_data)
override;
164 const unsigned integration_point)
const override
166 auto const& N =
_ip_data[integration_point].N;
169 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
174 std::vector<GlobalVector*>
const& ,
175 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
176 std::vector<double>& )
const override
184 std::vector<GlobalVector*>
const& ,
185 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
186 std::vector<double>& )
const override
194 std::vector<GlobalVector*>
const& ,
195 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
196 std::vector<double>& )
const override
204 std::vector<GlobalVector*>
const& ,
205 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
206 std::vector<double>& )
const override
214 std::vector<GlobalVector*>
const& ,
215 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
216 std::vector<double>& )
const override
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
unsigned getNumberOfPoints() const
virtual std::vector< double > const & getIntPtGasMolFracContaminant(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 & getIntPtLiquidMolFracContaminant(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
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 & getIntPtGasMolFracWater(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 > _gas_molar_fraction_water
typename ShapeMatricesType::GlobalDimNodalMatrixType GlobalDimNodalMatrixType
std::vector< double > _pressure_wetting
static const int temperature_matrix_index
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
std::vector< double > _saturation
NumLib::GenericIntegrationMethod const & _integration_method
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
typename LocalAssemblerTraits::LocalMatrix LocalMatrixType
typename LocalAssemblerTraits::LocalVector LocalVectorType
typename ShapeMatricesType::NodalMatrixType NodalMatrixType
std::vector< double > _liquid_molar_fraction_contaminant
static const int cap_pressure_size
typename ShapeMatricesType::NodalRowVectorType NodalRowVectorType
std::vector< double > const & getIntPtSaturation(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
ThermalTwoPhaseFlowWithPPLocalAssembler(MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, ThermalTwoPhaseFlowWithPPProcessData const &process_data)
std::vector< double > _gas_molar_fraction_contaminant
static const int nonwet_pressure_size
std::vector< double > const & getIntPtWettingPressure(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
static const int nonwet_pressure_matrix_index
std::vector< double > const & getIntPtGasMolFracWater(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
static const int temperature_size
static const int contaminant_matrix_index
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
static const int cap_pressure_matrix_index
std::vector< double > const & getIntPtLiquidMolFracContaminant(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
MeshLib::Element const & _element
static const int contaminant_size
typename ShapeMatricesType::NodalVectorType NodalVectorType
ThermalTwoPhaseFlowWithPPProcessData const & _process_data
std::vector< double > const & getIntPtGasMolFracContaminant(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
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)
const unsigned NUM_NODAL_DOF
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
NodalMatrixType const diffusion_operator
IntegrationPointData(NodalRowVectorType const &N_, GlobalDimNodalMatrixType const &dNdx_, double const &integration_weight_, NodalMatrixType const mass_operator_, NodalMatrixType const diffusion_operator_)
double const integration_weight
GlobalDimNodalMatrixType const dNdx
NodalMatrixType const mass_operator
NodalRowVectorType const N
EIGEN_MAKE_ALIGNED_OPERATOR_NEW
Vector< NNodes *NodalDOF > LocalVector
Matrix< NNodes *NodalDOF, NNodes *NodalDOF > LocalMatrix