29namespace TwoPhaseFlowWithPrho
31template <
typename NodalMatrixType>
69 std::vector<GlobalVector*>
const& x,
70 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
71 std::vector<double>& cache)
const = 0;
75 std::vector<GlobalVector*>
const& x,
76 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
77 std::vector<double>& cache)
const = 0;
80template <
typename ShapeFunction,
int GlobalDim>
102 bool const is_axially_symmetric,
116 unsigned const n_integration_points =
118 _ip_data.reserve(n_integration_points);
119 for (
unsigned ip = 0; ip < n_integration_points; ip++)
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;
137 void assemble(
double const t,
double const ,
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;
145 const unsigned integration_point)
const override
150 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
155 std::vector<GlobalVector*>
const& ,
156 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
157 std::vector<double>& )
const override
165 std::vector<GlobalVector*>
const& ,
166 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
167 std::vector<double>& )
const override
177 std::vector<ShapeMatrices, Eigen::aligned_allocator<ShapeMatrices>>
181 std::vector<IntegrationPointData<NodalMatrixType>,
182 Eigen::aligned_allocator<IntegrationPointData<NodalMatrixType>>>
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
unsigned getNumberOfPoints() 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
static const int cap_pressure_coeff_index
typename ShapeMatricesType::NodalVectorType NodalVectorType
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
static const int cap_pressure_matrix_index
static const int cap_pressure_size
static const int nonwet_pressure_matrix_index
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
static const int nonwet_pressure_size
TwoPhaseFlowWithPrhoLocalAssembler(MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, TwoPhaseFlowWithPrhoProcessData const &process_data)
NumLib::GenericIntegrationMethod const & _integration_method
std::vector< double > const & getIntPtNonWettingPressure(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
MeshLib::Element const & _element
static const int nonwet_pressure_coeff_index
std::vector< double > _saturation
std::vector< double > _pressure_nonwetting
used for secondary variable output
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
TwoPhaseFlowWithPrhoProcessData const & _process_data
typename LocalAssemblerTraits::LocalVector LocalVectorType
std::vector< double > const & getIntPtSaturation(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
typename LocalAssemblerTraits::LocalMatrix LocalMatrixType
typename ShapeMatricesType::NodalMatrixType NodalMatrixType
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
const unsigned NUM_NODAL_DOF
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
NodalMatrixType diffusionOperator
TwoPhaseFlowWithPrhoMaterialProperties & mat_property
double pressure_nonwetting
EIGEN_MAKE_ALIGNED_OPERATOR_NEW
IntegrationPointData(TwoPhaseFlowWithPrhoMaterialProperties &material_property_)
double integration_weight
NodalMatrixType massOperator
std::unique_ptr< TwoPhaseFlowWithPrhoMaterialProperties > _material
Vector< NNodes *NodalDOF > LocalVector
Matrix< NNodes *NodalDOF, NNodes *NodalDOF > LocalMatrix