33namespace WellboreSimulator
37template <
typename ShapeFunction,
int GlobalDim>
47 typename ShapeMatricesType::template VectorType<
NUM_NODAL_DOF *
48 ShapeFunction::NPOINTS>;
63 bool const is_axially_symmetric,
76 unsigned const n_integration_points =
78 _ip_data.reserve(n_integration_points);
80 auto const shape_matrices =
82 GlobalDim>(element, is_axially_symmetric,
89 auto const& liquid_phase = medium.
phase(
"AqueousLiquid");
91 for (
unsigned ip = 0; ip < n_integration_points; ip++)
94 shape_matrices[ip].N, shape_matrices[ip].dNdx,
96 shape_matrices[ip].integralMeasure *
97 shape_matrices[ip].detJ);
99 pos.setIntegrationPoint(ip);
116 .template value<double>(vars, pos, 0, 0);
121 .template value<double>(vars, pos, 0, 0);
123 _ip_data[ip].vapor_volume_fraction = 0;
124 _ip_data[ip].vapor_mass_flow_rate = 0;
125 _ip_data[ip].liquid_mass_flow_rate = 0;
130 void assemble(
double const t,
double const dt,
131 std::vector<double>
const& local_x,
132 std::vector<double>
const& local_x_prev,
133 std::vector<double>& local_M_data,
134 std::vector<double>& local_K_data,
135 std::vector<double>& local_b_data)
override;
145 double const liquid_water_density,
146 double const v_mix,
double const dryness,
147 double const C_0,
double const u_gu,
150 double const rho_mix =
151 alpha * vapor_water_density + (1 - alpha) * liquid_water_density;
154 dryness * liquid_water_density * rho_mix * v_mix -
155 alpha * C_0 * dryness * liquid_water_density * rho_mix * v_mix -
156 alpha * C_0 * (1 - dryness) * vapor_water_density * rho_mix *
158 alpha * vapor_water_density * liquid_water_density * u_gu;
162 double const liquid_water_density,
163 double const v_mix,
double const dryness,
164 double const C_0,
double const u_gu,
167 Jac(0) = dryness * liquid_water_density * v_mix *
168 (vapor_water_density - liquid_water_density) -
169 (C_0 * dryness * liquid_water_density +
170 C_0 * (1 - dryness) * vapor_water_density) *
171 (2 * alpha * vapor_water_density +
172 (1 - 2 * alpha) * liquid_water_density) *
174 vapor_water_density * liquid_water_density * u_gu;
178 const unsigned integration_point)
const override
180 auto const& N =
_ip_data[integration_point].N;
183 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
189 Eigen::VectorXd
const& ,
190 Eigen::VectorXd
const& )
override
192 auto const n_integration_points =
198 [](
double const s,
auto const& ip)
199 {
return s + ip.mix_density; }) /
200 n_integration_points;
204 Eigen::VectorXd
const& ,
205 double const ,
double const ,
208 unsigned const n_integration_points =
211 for (
unsigned ip = 0; ip < n_integration_points; ip++)
220 std::vector<GlobalVector*>
const& ,
221 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
222 std::vector<double>& cache)
const override
230 std::vector<GlobalVector*>
const& ,
231 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
232 std::vector<double>& cache)
const override
240 std::vector<GlobalVector*>
const& ,
241 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
242 std::vector<double>& cache)
const override
250 std::vector<GlobalVector*>
const& ,
251 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
252 std::vector<double>& cache)
const override
260 std::vector<GlobalVector*>
const& ,
261 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
262 std::vector<double>& cache)
const override
277 std::vector<IpData, Eigen::aligned_allocator<IpData>>
_ip_data;
Medium * getMedium(std::size_t element_id)
Phase const & phase(std::size_t index) const
double liquid_phase_pressure
Eigen::Vector3d const & asEigenVector3d() const
virtual const Node * getNode(unsigned idx) const =0
std::size_t getID() const
Returns the ID of the element.
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
unsigned getNumberOfPoints() const
void setElementID(std::size_t element_id)
virtual std::vector< double > const & getIntPtVaporVolumeFraction(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
MeshLib::Element const & _element
WellboreSimulatorProcessData const & _process_data
virtual std::vector< double > const & getIntPtTemperature(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
typename ShapeMatricesType::NodalVectorType NodalVectorType
WellboreSimulatorFEM(MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, WellboreSimulatorProcessData const &process_data)
typename ShapeMatricesType::template MatrixType< NUM_NODAL_DOF *ShapeFunction::NPOINTS, NUM_NODAL_DOF *ShapeFunction::NPOINTS > LocalMatrixType
static const int velocity_index
Eigen::Matrix< double, jacobian_residual_size, 1 > UnknownVector
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
static const int velocity_size
typename ShapeMatricesType::GlobalDimNodalMatrixType GlobalDimNodalMatrixType
static const int pressure_index
void calculateJacobian(double const alpha, double const vapor_water_density, double const liquid_water_density, double const v_mix, double const dryness, double const C_0, double const u_gu, JacobianMatrix &Jac)
bool const _is_axially_symmetric
Eigen::Matrix< double, jacobian_residual_size, 1 > ResidualVector
NumLib::GenericIntegrationMethod const & _integration_method
static const int enthalpy_index
Eigen::Matrix< double, jacobian_residual_size, jacobian_residual_size, Eigen::RowMajor > JacobianMatrix
Eigen::Vector3d _element_direction
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
static const int enthalpy_size
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
typename ShapeMatricesType::NodalRowVectorType NodalRowVectorType
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
virtual std::vector< double > const & getIntPtLiquidMassFlowRate(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
typename ShapeMatricesType::template VectorType< NUM_NODAL_DOF * ShapeFunction::NPOINTS > LocalVectorType
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
void postTimestepConcrete(Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override
void computeSecondaryVariableConcrete(double const, double const, Eigen::VectorXd const &, Eigen::VectorXd const &) override
void calculateResidual(double const alpha, double const vapor_water_density, double const liquid_water_density, double const v_mix, double const dryness, double const C_0, double const u_gu, ResidualVector &res)
static const int pressure_size
virtual std::vector< double > const & getIntPtDryness(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
virtual std::vector< double > const & getIntPtVaporMassFlowRate(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
static int const jacobian_residual_size
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_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)
const unsigned NUM_NODAL_DOF
std::vector< double > const & getIntegrationPointScalarData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)
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
virtual Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > getNodalValuesOnElement(MeshLib::Element const &element, double const t) const
Returns a matrix of values for all nodes of the given element.
double liquid_mass_flow_rate
double vapor_volume_fraction
double vapor_mass_flow_rate
MaterialPropertyLib::MaterialSpatialDistributionMap media_map
ParameterLib::Parameter< double > const & well_ref_enthalpy
MeshLib::PropertyVector< double > * mesh_prop_density
ParameterLib::Parameter< double > const & well_ref_pressure