28namespace TwoPhaseFlowWithPP
30template <
typename NodalRowVectorType,
typename GlobalDimNodalMatrixType,
31 typename NodalMatrixType>
35 GlobalDimNodalMatrixType dNdx_,
36 double const& integration_weight_)
38 dNdx(std::move(dNdx_)),
43 NodalRowVectorType
const N;
44 GlobalDimNodalMatrixType
const dNdx;
58 std::vector<GlobalVector*>
const& x,
59 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
60 std::vector<double>& cache)
const = 0;
64 std::vector<GlobalVector*>
const& x,
65 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
66 std::vector<double>& cache)
const = 0;
69template <
typename ShapeFunction,
int GlobalDim>
94 bool const is_axially_symmetric,
104 unsigned const n_integration_points =
106 _ip_data.reserve(n_integration_points);
107 auto const shape_matrices =
109 GlobalDim>(element, is_axially_symmetric,
111 for (
unsigned ip = 0; ip < n_integration_points; ip++)
113 auto const& sm = shape_matrices[ip];
116 sm.integralMeasure * sm.detJ *
123 int const process_id)
override;
125 void assemble(
double const t,
double const dt,
126 std::vector<double>
const& local_x,
127 std::vector<double>
const& local_x_prev,
128 std::vector<double>& local_M_data,
129 std::vector<double>& local_K_data,
130 std::vector<double>& local_b_data)
override;
133 const unsigned integration_point)
const override
135 auto const& N =
_ip_data[integration_point].N;
138 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
143 std::vector<GlobalVector*>
const& ,
144 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
145 std::vector<double>& )
const override
153 std::vector<GlobalVector*>
const& ,
154 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
155 std::vector<double>& )
const override
constexpr double getWeight() const
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 & 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
static const int cap_pressure_coeff_index
typename ShapeMatricesType::NodalRowVectorType NodalRowVectorType
TwoPhaseFlowWithPPLocalAssembler(MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, TwoPhaseFlowWithPPProcessData const &process_data)
typename LocalAssemblerTraits::LocalMatrix LocalMatrixType
typename LocalAssemblerTraits::LocalVector LocalVectorType
static const int nonwet_pressure_matrix_index
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType > > > _ip_data
static const int nonwet_pressure_coeff_index
static const int cap_pressure_matrix_index
static const int nonwet_pressure_size
typename ShapeMatricesType::NodalMatrixType NodalMatrixType
typename ShapeMatricesType::GlobalDimNodalMatrixType GlobalDimNodalMatrixType
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
MeshLib::Element const & _element
std::vector< double > const & getIntPtWetPressure(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
std::vector< double > _pressure_wet
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
TwoPhaseFlowWithPPProcessData const & _process_data
static const int cap_pressure_size
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
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
typename ShapeMatricesType::NodalVectorType NodalVectorType
std::vector< double > _saturation
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
IntegrationPointData(NodalRowVectorType N_, GlobalDimNodalMatrixType dNdx_, double const &integration_weight_)
GlobalDimNodalMatrixType const dNdx
NodalRowVectorType const N
const double integration_weight
EIGEN_MAKE_ALIGNED_OPERATOR_NEW
Vector< NNodes *NodalDOF > LocalVector
Matrix< NNodes *NodalDOF, NNodes *NodalDOF > LocalMatrix