OGS
WellboreSimulatorFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <Eigen/Dense>
14#include <numeric>
15#include <vector>
16
30
31namespace ProcessLib
32{
33namespace WellboreSimulator
34{
35const unsigned NUM_NODAL_DOF = 3;
36
37template <typename ShapeFunction, int GlobalDim>
39{
42
43 using LocalMatrixType = typename ShapeMatricesType::template MatrixType<
44 NUM_NODAL_DOF * ShapeFunction::NPOINTS,
45 NUM_NODAL_DOF * ShapeFunction::NPOINTS>;
47 typename ShapeMatricesType::template VectorType<NUM_NODAL_DOF *
48 ShapeFunction::NPOINTS>;
49
52
57
58public:
60 MeshLib::Element const& element,
61 std::size_t const /*local_matrix_size*/,
62 NumLib::GenericIntegrationMethod const& integration_method,
63 bool const is_axially_symmetric,
64 WellboreSimulatorProcessData const& process_data)
65 : _element(element),
66 _integration_method(integration_method),
67 _is_axially_symmetric(is_axially_symmetric),
68 _process_data(process_data)
69 {
70 // calculate the element direction vector
71 auto const& p0 = element.getNode(0)->asEigenVector3d();
72 auto const& p1 = element.getNode(1)->asEigenVector3d();
73
74 _element_direction = (p1 - p0).normalized();
75
76 unsigned const n_integration_points =
78 _ip_data.reserve(n_integration_points);
79
80 auto const shape_matrices =
82 GlobalDim>(element, is_axially_symmetric,
84
87 auto const& medium =
89 auto const& liquid_phase = medium.phase("AqueousLiquid");
90
91 for (unsigned ip = 0; ip < n_integration_points; ip++)
92 {
93 _ip_data.emplace_back(
94 shape_matrices[ip].N, shape_matrices[ip].dNdx,
96 shape_matrices[ip].integralMeasure *
97 shape_matrices[ip].detJ);
98
99 pos.setIntegrationPoint(ip);
101
102 NodalVectorType ref_p =
104 _element, 0);
105 NodalVectorType ref_h =
107 _element, 0);
108
109 vars.liquid_phase_pressure = _ip_data[ip].N.dot(ref_p);
110 vars.enthalpy = _ip_data[ip].N.dot(ref_h);
111
112 // .initialValue
113 _ip_data[ip].temperature =
114 liquid_phase
116 .template value<double>(vars, pos, 0, 0);
117 vars.temperature = _ip_data[ip].temperature;
118 _ip_data[ip].mix_density =
119 liquid_phase
121 .template value<double>(vars, pos, 0, 0);
122 _ip_data[ip].dryness = 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;
126 _ip_data[ip].pushBackState();
127 }
128 }
129
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;
136
137 static int const jacobian_residual_size = 1;
138 using ResidualVector = Eigen::Matrix<double, jacobian_residual_size, 1>;
140 Eigen::Matrix<double, jacobian_residual_size, jacobian_residual_size,
141 Eigen::RowMajor>;
142 using UnknownVector = Eigen::Matrix<double, jacobian_residual_size, 1>;
143
144 void calculateResidual(double const alpha, double const vapor_water_density,
145 double const liquid_water_density,
146 double const v_mix, double const dryness,
147 double const C_0, double const u_gu,
148 ResidualVector& res)
149 {
150 double const rho_mix =
151 alpha * vapor_water_density + (1 - alpha) * liquid_water_density;
152
153 res(0) =
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 *
157 v_mix -
158 alpha * vapor_water_density * liquid_water_density * u_gu;
159 }
160
161 void calculateJacobian(double const alpha, double const vapor_water_density,
162 double const liquid_water_density,
163 double const v_mix, double const dryness,
164 double const C_0, double const u_gu,
165 JacobianMatrix& Jac)
166 {
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) *
173 v_mix -
174 vapor_water_density * liquid_water_density * u_gu;
175 }
176
177 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
178 const unsigned integration_point) const override
179 {
180 auto const& N = _ip_data[integration_point].N;
181
182 // assumes N is stored contiguously in memory
183 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
184 }
185
187 double const /*t*/,
188 double const /*dt*/,
189 Eigen::VectorXd const& /*local_x*/,
190 Eigen::VectorXd const& /*local_x_dot*/) override
191 {
192 auto const n_integration_points =
194 auto const ele_id = _element.getID();
195
197 std::accumulate(_ip_data.begin(), _ip_data.end(), 0.,
198 [](double const s, auto const& ip)
199 { return s + ip.mix_density; }) /
200 n_integration_points;
201 }
202
203 void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
204 Eigen::VectorXd const& /*local_x_dot*/,
205 double const /*t*/, double const /*dt*/,
206 int const /*process_id*/) override
207 {
208 unsigned const n_integration_points =
210
211 for (unsigned ip = 0; ip < n_integration_points; ip++)
212 {
213 _ip_data[ip].pushBackState();
214 }
215 }
216
217protected:
218 virtual std::vector<double> const& getIntPtVaporMassFlowRate(
219 const double /*t*/,
220 std::vector<GlobalVector*> const& /*x*/,
221 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
222 std::vector<double>& cache) const override
223 {
226 }
227
228 virtual std::vector<double> const& getIntPtLiquidMassFlowRate(
229 const double /*t*/,
230 std::vector<GlobalVector*> const& /*x*/,
231 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
232 std::vector<double>& cache) const override
233 {
236 }
237
238 virtual std::vector<double> const& getIntPtTemperature(
239 const double /*t*/,
240 std::vector<GlobalVector*> const& /*x*/,
241 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
242 std::vector<double>& cache) const override
243 {
246 }
247
248 virtual std::vector<double> const& getIntPtDryness(
249 const double /*t*/,
250 std::vector<GlobalVector*> const& /*x*/,
251 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
252 std::vector<double>& cache) const override
253 {
255 _ip_data, &IpData::dryness, cache);
256 }
257
258 virtual std::vector<double> const& getIntPtVaporVolumeFraction(
259 const double /*t*/,
260 std::vector<GlobalVector*> const& /*x*/,
261 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
262 std::vector<double>& cache) const override
263 {
266 }
267
272
273 Eigen::Vector3d _element_direction;
274
275 using IpData =
277 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
278
279protected:
280 static const int pressure_index = 0;
281 static const int velocity_index = ShapeFunction::NPOINTS;
282 static const int enthalpy_index = 2 * ShapeFunction::NPOINTS;
283
284 static const int pressure_size = ShapeFunction::NPOINTS;
285 static const int velocity_size = ShapeFunction::NPOINTS;
286 static const int enthalpy_size = ShapeFunction::NPOINTS;
287};
288
289} // namespace WellboreSimulator
290} // namespace ProcessLib
291
Phase const & phase(std::size_t index) const
Definition Medium.cpp:33
Eigen::Vector3d const & asEigenVector3d() const
Definition Point3d.h:64
double getWeight() const
virtual const Node * getNode(unsigned idx) const =0
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) 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
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
Eigen::Matrix< double, jacobian_residual_size, 1 > UnknownVector
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
typename ShapeMatricesType::GlobalDimNodalMatrixType GlobalDimNodalMatrixType
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)
Eigen::Matrix< double, jacobian_residual_size, 1 > ResidualVector
NumLib::GenericIntegrationMethod const & _integration_method
Eigen::Matrix< double, jacobian_residual_size, jacobian_residual_size, Eigen::RowMajor > JacobianMatrix
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
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)
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
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)
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.
Definition Parameter.h:164
MaterialPropertyLib::MaterialSpatialDistributionMap media_map