OGS
ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, int GlobalDim>
class ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >

Definition at line 31 of file WellboreSimulatorFEM.h.

#include <WellboreSimulatorFEM.h>

Inheritance diagram for ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >:
[legend]

Public Types

using ResidualVector = Eigen::Matrix<double, jacobian_residual_size, 1>
using JacobianMatrix
using UnknownVector = Eigen::Matrix<double, jacobian_residual_size, 1>

Public Member Functions

 WellboreSimulatorFEM (MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, WellboreSimulatorProcessData const &process_data)
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 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)
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::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
void computeSecondaryVariableConcrete (double const, double const, Eigen::VectorXd const &, Eigen::VectorXd const &) override
void postTimestepConcrete (Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override
Public Member Functions inherited from ProcessLib::WellboreSimulator::WellboreSimulatorLocalAssemblerInterface
 WellboreSimulatorLocalAssemblerInterface ()=default
Public Member Functions inherited from ProcessLib::LocalAssemblerInterface
virtual ~LocalAssemblerInterface ()=default
virtual void setInitialConditions (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, int const process_id)
virtual void initialize (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
virtual void preAssemble (double const, double const, std::vector< double > const &)
virtual void assembleForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
virtual void assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
virtual void assembleWithJacobianForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
virtual void computeSecondaryVariable (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
virtual void preTimestep (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
virtual void postTimestep (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
void postNonLinearSolver (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< double > const &) const
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double > > const &) const
 Fits to staggered scheme.
virtual int getNumberOfVectorElementsForDeformation () const
Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default

Static Public Attributes

static int const jacobian_residual_size = 1

Protected Types

using IpData

Protected Member Functions

virtual std::vector< double > const & getIntPtVaporMassFlowRate (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
virtual std::vector< double > const & getIntPtLiquidMassFlowRate (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
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 & getIntPtVaporVolumeFraction (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override

Protected Attributes

MeshLib::Element const & _element
NumLib::GenericIntegrationMethod const & _integration_method
bool const _is_axially_symmetric
WellboreSimulatorProcessData const & _process_data
Eigen::Vector3d _element_direction
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data

Static Protected Attributes

static const int pressure_index = 0
static const int velocity_index = ShapeFunction::NPOINTS
static const int enthalpy_index = 2 * ShapeFunction::NPOINTS
static const int pressure_size = ShapeFunction::NPOINTS
static const int velocity_size = ShapeFunction::NPOINTS
static const int enthalpy_size = ShapeFunction::NPOINTS

Private Types

using ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
using LocalMatrixType
using LocalVectorType
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
using NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
using GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
using GlobalDimNodalMatrixType
using GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType

Member Typedef Documentation

◆ GlobalDimMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType
private

Definition at line 49 of file WellboreSimulatorFEM.h.

◆ GlobalDimNodalMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::GlobalDimNodalMatrixType
private
Initial value:
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType

Definition at line 47 of file WellboreSimulatorFEM.h.

◆ GlobalDimVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
private

Definition at line 46 of file WellboreSimulatorFEM.h.

◆ IpData

◆ JacobianMatrix

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::JacobianMatrix
Initial value:

Definition at line 131 of file WellboreSimulatorFEM.h.

◆ LocalMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::LocalMatrixType
private
Initial value:
typename ShapeMatricesType::template MatrixType<
NUM_NODAL_DOF * ShapeFunction::NPOINTS,
NUM_NODAL_DOF * ShapeFunction::NPOINTS>
const unsigned NUM_NODAL_DOF

Definition at line 36 of file WellboreSimulatorFEM.h.

◆ LocalVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::LocalVectorType
private
Initial value:
typename ShapeMatricesType::template VectorType<NUM_NODAL_DOF *
ShapeFunction::NPOINTS>

Definition at line 39 of file WellboreSimulatorFEM.h.

◆ NodalRowVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
private

Definition at line 44 of file WellboreSimulatorFEM.h.

◆ NodalVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::NodalVectorType = typename ShapeMatricesType::NodalVectorType
private

Definition at line 43 of file WellboreSimulatorFEM.h.

◆ ResidualVector

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::ResidualVector = Eigen::Matrix<double, jacobian_residual_size, 1>

Definition at line 130 of file WellboreSimulatorFEM.h.

◆ ShapeMatrices

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
private

Definition at line 34 of file WellboreSimulatorFEM.h.

◆ ShapeMatricesType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
private

Definition at line 33 of file WellboreSimulatorFEM.h.

◆ UnknownVector

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::UnknownVector = Eigen::Matrix<double, jacobian_residual_size, 1>

Definition at line 134 of file WellboreSimulatorFEM.h.

Constructor & Destructor Documentation

◆ WellboreSimulatorFEM()

template<typename ShapeFunction, int GlobalDim>
ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::WellboreSimulatorFEM ( MeshLib::Element const & element,
std::size_t const ,
NumLib::GenericIntegrationMethod const & integration_method,
bool const is_axially_symmetric,
WellboreSimulatorProcessData const & process_data )
inline

Definition at line 52 of file WellboreSimulatorFEM.h.

62 {
63 // calculate the element direction vector
64 auto const& p0 = element.getNode(0)->asEigenVector3d();
65 auto const& p1 = element.getNode(1)->asEigenVector3d();
66
68
69 unsigned const n_integration_points =
70 _integration_method.getNumberOfPoints();
72
73 auto const shape_matrices =
77
79 pos.setElementID(_element.getID());
80 auto const& medium =
81 *_process_data.media_map.getMedium(_element.getID());
82 auto const& liquid_phase = medium.phase("AqueousLiquid");
83
84 for (unsigned ip = 0; ip < n_integration_points; ip++)
85 {
86 _ip_data.emplace_back(
88 _integration_method.getWeightedPoint(ip).getWeight() *
89 shape_matrices[ip].integralMeasure *
90 shape_matrices[ip].detJ);
91
93
95 _process_data.well_ref_pressure.getNodalValuesOnElement(
96 _element, 0);
98 _process_data.well_ref_enthalpy.getNodalValuesOnElement(
99 _element, 0);
100
101 vars.liquid_phase_pressure = _ip_data[ip].N.dot(ref_p);
102 vars.enthalpy = _ip_data[ip].N.dot(ref_h);
103
104 // .initialValue
105 _ip_data[ip].temperature =
108 .template value<double>(vars, pos, 0, 0);
109 vars.temperature = _ip_data[ip].temperature;
110 _ip_data[ip].mix_density =
113 .template value<double>(vars, pos, 0, 0);
114 _ip_data[ip].dryness = 0;
115 _ip_data[ip].vapor_volume_fraction = 0;
116 _ip_data[ip].vapor_mass_flow_rate = 0;
117 _ip_data[ip].liquid_mass_flow_rate = 0;
118 _ip_data[ip].pushBackState();
119 }
120 }
typename ShapeMatricesType::NodalVectorType NodalVectorType
NumLib::GenericIntegrationMethod const & _integration_method
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data

References _element, _element_direction, _integration_method, _ip_data, _is_axially_symmetric, _process_data, MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::enthalpy, MeshLib::Element::getNode(), NumLib::initShapeMatrices(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, ParameterLib::SpatialPosition::setElementID(), MaterialPropertyLib::temperature, and MaterialPropertyLib::VariableArray::temperature.

Member Function Documentation

◆ assemble()

template<typename ShapeFunction, int GlobalDim>
void ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::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 )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 21 of file WellboreSimulatorFEM-impl.h.

25{
26 auto const local_matrix_size = local_x.size();
27
29
36
37 // Get block matrices
40
45
48
53
56
60
61 unsigned const n_integration_points =
62 _integration_method.getNumberOfPoints();
63
65 pos.setElementID(_element.getID());
66
67 auto const& b = _process_data.specific_body_force;
68
70
71 // Get material properties
72 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
73 auto const& liquid_phase = medium.phase("AqueousLiquid");
74 auto const& gas_phase = medium.phase("Gas");
75
76 // Get wellbore parameters
77 // casing thickness
78 auto const t_ca = _process_data.wellbore.casing_thickness(t, pos)[0];
79 // wellbore radius
80 auto const r_w = _process_data.wellbore.diameter(t, pos)[0] / 2;
81
82 // pipe thickness
83 auto const t_p = _process_data.wellbore.pipe_thickness(t, pos)[0];
84
85 // roughness of the wellbore
86 auto const xi = _process_data.wellbore.roughness(t, pos)[0];
87 // pipe outer radius
88 auto const r_o = r_w - t_ca;
89 // pipe inner radius
90 auto const r_i = r_o - t_p;
91
92 // get reservoir properties
94 _process_data.reservoir_properties.temperature.getNodalValuesOnElement(
95 _element, t);
97 _process_data.reservoir_properties.pressure.getNodalValuesOnElement(
98 _element, t);
100 _process_data.productivity_index.getNodalValuesOnElement(_element, t);
101 auto const k_r =
102 _process_data.reservoir_properties.thermal_conductivity(t, pos)[0];
103 auto const rho_r = _process_data.reservoir_properties.density(t, pos)[0];
104 auto const c_r =
105 _process_data.reservoir_properties.specific_heat_capacity(t, pos)[0];
106
107 for (unsigned ip(0); ip < n_integration_points; ip++)
108 {
109 auto& ip_data = _ip_data[ip];
110 auto const& N = ip_data.N;
111 auto const& dNdx = ip_data.dNdx;
112 auto const& w = ip_data.integration_weight;
113 auto& mix_density = ip_data.mix_density;
114 auto& temperature = ip_data.temperature;
115 auto& steam_mass_frac = ip_data.dryness;
116 auto& vapor_volume_frac = ip_data.vapor_volume_fraction;
117 auto& vapor_mass_flowrate = ip_data.vapor_mass_flow_rate;
118 auto& liquid_mass_flowrate = ip_data.liquid_mass_flow_rate;
119
120 pos = {
121 std::nullopt, _element.getID(),
124 _element, N))};
125
126 double p_int_pt = 0.0;
127 double v_int_pt = 0.0;
128 double h_int_pt = 0.0;
129
131 h_int_pt);
132
133 double p_prev_int_pt = 0.0;
134 double v_prev_int_pt = 0.0;
135 double h_prev_int_pt = 0.0;
136
139
140 double vdot_int_pt = (v_int_pt - v_prev_int_pt) / dt;
141
142 // calculate fluid properties
143
144 const double pi = std::numbers::pi;
145
146 vars.liquid_phase_pressure = p_int_pt;
147 vars.enthalpy = h_int_pt;
148
149 double liquid_water_density =
152 .template value<double>(vars, pos, t, dt);
153 double const vapour_water_density =
156 .template value<double>(vars, pos, t, dt);
157
158 double const h_sat_liq_w =
160 .property(
162 .template value<double>(vars, pos, t, dt);
163 double const h_sat_vap_w =
165 .property(
167 .template value<double>(vars, pos, t, dt);
168
169 // TODO add a function to calculate dryness with constrain of 0
170 // to 1.
171 double const dryness = std::max(
174
175 double const T_int_pt =
176 (dryness == 0)
179 .template value<double>(vars, pos, t, dt)
180 : gas_phase
185 vars.temperature = T_int_pt;
186
187 // For the calculation of the void fraction of vapour,
188 // see Rohuani, Z., and E. Axelsson. "Calculation of volume void
189 // fraction in a subcooled and quality region." International
190 // Journal of Heat and Mass Transfer 17 (1970): 383-393.
191
192 // profile parameter of drift flux
193 double C_0 = 1 + 0.12 * (1 - dryness);
194
195 // For the surface tension calculation, see
196 // Cooper, J. R., and R. B. Dooley. "IAPWS release on surface
197 // tension of ordinary water substance." International Association
198 // for the Properties of Water and Steam (1994).
199 double const sigma_gl = 0.2358 *
200 std::pow((1 - T_int_pt / 647.096), 1.256) *
201 (1 - 0.625 * (1 - T_int_pt / 647.096));
202 // drift flux velocity
203 double const u_gu =
204 1.18 * (1 - dryness) *
205 std::pow((9.81) * sigma_gl *
207 0.25) /
209
210 // solving void fraction of vapor: Rouhani-Axelsson
211 double alpha = 0;
212 if (dryness != 0)
213 {
214 // Local Newton solver
215 using LocalJacobianMatrix =
220
222
224 {
227 u_gu, residual);
228 };
229
231 {
234 dryness, C_0, u_gu,
235 jacobian); // for solution dependent Jacobians
236 };
237
238 auto const update_solution =
239 [&](LocalUnknownVector const& increment)
240 {
241 // increment solution vectors
242 alpha += increment[0];
243 };
244
245 const int maximum_iterations(20);
246 const double residuum_tolerance(1.e-10);
247 const double increment_tolerance(0);
248
253
254 auto const success_iterations = newton_solver.solve(J_loc);
255
257 {
258 WARN(
259 "Attention! Steam void fraction has not been correctly "
260 "calculated!");
261 }
262 }
263
265
266 if (alpha == 0)
267 {
271 .template value<double>(vars, pos, t, dt);
272 }
273
276
277 auto& mix_density_prev = ip_data.mix_density_prev;
278 vars.density = mix_density;
279
280 auto const rho_dot = (mix_density - mix_density_prev) / dt;
281
282 double const liquid_water_velocity_act =
283 (alpha == 0) ? v_int_pt
284 : (1 - dryness) * mix_density * v_int_pt /
286 double const vapor_water_velocity_act =
287 (alpha == 0) ? 0
290
292 pi * r_i * r_i * alpha;
293
296 (1 - alpha);
297
298 // Slip parameter between two phases,
299 // see Akbar, Somaieh, N. Fathianpour, and Rafid Al Khoury. "A finite
300 // element model for high enthalpy two-phase flow in geothermal
301 // wellbores." Renewable Energy 94 (2016): 223-236.
302 double const gamma =
304 (1 - alpha) /
306 (1 - alpha * C_0) * liquid_water_density),
307 2) *
308 std::pow((C_0 - 1) * v_int_pt + u_gu, 2);
309
310 double const miu =
312 .template value<double>(vars, pos, t, dt);
313 double const Re = mix_density * v_int_pt * 2 * r_i / miu;
314
315 // Wall friction coefficient,
316 // Musivand Arzanfudi, Mehdi, and Rafid Al‐Khoury. "A compressible
317 // two‐fluid multiphase model for CO2 leakage through a wellbore."
318 // International Journal for Numerical Methods in Fluids 77.8 (2015):
319 // 477-507.
320 double f = 0.0;
321 if (Re > 10 && Re <= 2400)
322 f = 16 / Re;
323 else if (Re > 2400)
324 f = std::pow(std::log(xi / 3.7 / r_i) -
325 5.02 / Re * std::log(xi / 3.7 / r_i + 13 / Re),
326 -2) /
327 16;
328
329 double Q_hx = 0;
330 double const T_r_int_pt = N.dot(T_r);
331 // conductive heat exchange between wellbore and formation
332 if (_process_data.has_heat_exchange_with_formation)
333 {
334 // See Zhang, Pan, Pruess, Finsterle (2011). A time-convolution
335 // approach for modeling heat exchange between a wellbore and
336 // surrounding formation. Geothermics 40, 261-266.
337 const double alpha_r = k_r / rho_r / c_r;
338 const double t_d = alpha_r * t / (r_i * r_i);
339
340 double beta;
341 if (t_d < 2.8)
342 beta = std::pow((pi * t_d), -0.5) + 0.5 -
343 0.25 * std::pow((t_d / pi), 0.5) + 0.125 * t_d;
344 else
345 beta = 2 * (1 / (std::log(4 * t_d) - 2 * 0.57722) -
346 0.57722 /
347 std::pow((std::log(4 * t_d) - 2 * 0.57722), 2));
348
349 const double P_c = 2 * pi * r_i;
350 Q_hx = P_c * k_r * (T_r_int_pt - T_int_pt) / r_i * beta;
351 }
352
353 // mass exchange with reservoir
354 double const p_r_int_pt = N.dot(p_r);
355 double const PI_int_pt = N.dot(PI);
356 double Q_mx = PI_int_pt * (p_int_pt - p_r_int_pt);
357
358 // advective momentum and energy exchange with reservoir due to the mass
359 // exchange
360 double Q_mom = 0;
361 double Q_ene = 0;
362 if (Q_mx != 0)
363 {
364 Q_mom = Q_mx * v_int_pt;
365 // Only single-phase liquid condition is considered now
366 // TODO: update the two-phase flowing enthalpy from the feed zone.
367 vars.liquid_phase_pressure = p_r_int_pt;
368 vars.temperature = T_r_int_pt;
369 double const h_fres =
372 .template value<double>(vars, pos, t, dt);
373 Q_ene = Q_mx * h_fres;
374 }
375
376 // M matrix assembly
377 Mvv.noalias() += w * N.transpose() * mix_density * N;
378
379 Mhp.noalias() += -w * N.transpose() * N;
380 Mhh.noalias() += w * N.transpose() * mix_density * N;
381
382 // K matrix assembly
383 Kpv.noalias() += w * dNdx.transpose() * N * mix_density;
384
385 Kvp.noalias() += w * N.transpose() * dNdx;
386 Kvv.noalias() += w * N.transpose() * rho_dot * N;
387
388 Khh.noalias() += w * N.transpose() * mix_density * v_int_pt * dNdx;
389
390 // b matrix assembly
391 Bp.noalias() += w * N.transpose() * rho_dot + w * N.transpose() * Q_mx;
392
393 Bv.noalias() +=
394 w * dNdx.transpose() * mix_density * v_int_pt * v_int_pt +
395 w * dNdx.transpose() * gamma -
396 w * N.transpose() * f * mix_density * std::abs(v_int_pt) *
397 v_int_pt / (4 * r_i) -
398 w * N.transpose() * Q_mom;
399
400 Bh.noalias() +=
401 -1 / 2 * w * N.transpose() * rho_dot * v_int_pt * v_int_pt -
402 w * N.transpose() * mix_density * v_int_pt * vdot_int_pt +
403 1 / 2 * w * dNdx.transpose() * mix_density * v_int_pt * v_int_pt *
404 v_int_pt +
405 w * N.transpose() * (Q_hx / pi / r_i / r_i) -
406 w * N.transpose() * Q_ene;
407
408 if (_process_data.has_gravity)
409 {
411 N.transpose() * b * w * _element_direction[2];
412
413 Bv.noalias() += gravity_operator * mix_density;
414 Bh.noalias() += gravity_operator * mix_density * v_int_pt;
415 }
416 }
417
418 // debugging
419 // std::string sep = "\n----------------------------------------\n";
420 // Eigen::IOFormat CleanFmt(6, 0, ", ", "\n", "[", "]");
421 // std::cout << local_M.format(CleanFmt) << sep;
422 // std::cout << local_K.format(CleanFmt) << sep;
423 // std::cout << local_b.format(CleanFmt) << sep;
424}
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:34
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)
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)
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)

References _element, _element_direction, _integration_method, _ip_data, _process_data, calculateJacobian(), calculateResidual(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, MaterialPropertyLib::enthalpy, MaterialPropertyLib::VariableArray::enthalpy, enthalpy_index, NumLib::interpolateCoordinates(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, ProcessLib::WellboreSimulator::NUM_NODAL_DOF, PI, pressure_index, MaterialPropertyLib::saturation_density, MaterialPropertyLib::saturation_enthalpy, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), NumLib::NewtonRaphson< LinearSolver, JacobianMatrixUpdate, ResidualUpdate, SolutionUpdate >::solve(), MaterialPropertyLib::temperature, MaterialPropertyLib::VariableArray::temperature, velocity_index, MaterialPropertyLib::viscosity, and WARN().

◆ calculateJacobian()

template<typename ShapeFunction, int GlobalDim>
void ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::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 )
inline

Definition at line 153 of file WellboreSimulatorFEM.h.

Referenced by assemble().

◆ calculateResidual()

template<typename ShapeFunction, int GlobalDim>
void ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::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 )
inline

Definition at line 136 of file WellboreSimulatorFEM.h.

Referenced by assemble().

◆ computeSecondaryVariableConcrete()

template<typename ShapeFunction, int GlobalDim>
void ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::computeSecondaryVariableConcrete ( double const ,
double const ,
Eigen::VectorXd const & ,
Eigen::VectorXd const &  )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 178 of file WellboreSimulatorFEM.h.

183 {
184 auto const n_integration_points =
185 _integration_method.getNumberOfPoints();
186 auto const ele_id = _element.getID();
187
188 (*_process_data.mesh_prop_density)[ele_id] =
189 std::accumulate(_ip_data.begin(), _ip_data.end(), 0.,
190 [](double const s, auto const& ip)
191 { return s + ip.mix_density; }) /
193 }

References _element, _integration_method, _ip_data, and _process_data.

◆ getIntPtDryness()

template<typename ShapeFunction, int GlobalDim>
virtual std::vector< double > const & ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::getIntPtDryness ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverrideprotectedvirtual

◆ getIntPtLiquidMassFlowRate()

template<typename ShapeFunction, int GlobalDim>
virtual std::vector< double > const & ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::getIntPtLiquidMassFlowRate ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverrideprotectedvirtual

◆ getIntPtTemperature()

template<typename ShapeFunction, int GlobalDim>
virtual std::vector< double > const & ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::getIntPtTemperature ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverrideprotectedvirtual

◆ getIntPtVaporMassFlowRate()

template<typename ShapeFunction, int GlobalDim>
virtual std::vector< double > const & ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::getIntPtVaporMassFlowRate ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverrideprotectedvirtual

◆ getIntPtVaporVolumeFraction()

template<typename ShapeFunction, int GlobalDim>
virtual std::vector< double > const & ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::getIntPtVaporVolumeFraction ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverrideprotectedvirtual

◆ getShapeMatrix()

template<typename ShapeFunction, int GlobalDim>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::getShapeMatrix ( const unsigned integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 169 of file WellboreSimulatorFEM.h.

171 {
172 auto const& N = _ip_data[integration_point].N;
173
174 // assumes N is stored contiguously in memory
175 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
176 }

References _ip_data.

◆ postTimestepConcrete()

template<typename ShapeFunction, int GlobalDim>
void ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::postTimestepConcrete ( Eigen::VectorXd const & ,
Eigen::VectorXd const & ,
double const ,
double const ,
int const  )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 195 of file WellboreSimulatorFEM.h.

199 {
200 unsigned const n_integration_points =
201 _integration_method.getNumberOfPoints();
202
203 for (unsigned ip = 0; ip < n_integration_points; ip++)
204 {
205 _ip_data[ip].pushBackState();
206 }
207 }

References _integration_method, and _ip_data.

Member Data Documentation

◆ _element

template<typename ShapeFunction, int GlobalDim>
MeshLib::Element const& ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::_element
protected

◆ _element_direction

template<typename ShapeFunction, int GlobalDim>
Eigen::Vector3d ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::_element_direction
protected

Definition at line 265 of file WellboreSimulatorFEM.h.

Referenced by WellboreSimulatorFEM(), and assemble().

◆ _integration_method

◆ _ip_data

◆ _is_axially_symmetric

template<typename ShapeFunction, int GlobalDim>
bool const ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::_is_axially_symmetric
protected

Definition at line 262 of file WellboreSimulatorFEM.h.

Referenced by WellboreSimulatorFEM().

◆ _process_data

template<typename ShapeFunction, int GlobalDim>
WellboreSimulatorProcessData const& ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::_process_data
protected

◆ enthalpy_index

template<typename ShapeFunction, int GlobalDim>
const int ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::enthalpy_index = 2 * ShapeFunction::NPOINTS
staticprotected

Definition at line 274 of file WellboreSimulatorFEM.h.

Referenced by assemble().

◆ enthalpy_size

template<typename ShapeFunction, int GlobalDim>
const int ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::enthalpy_size = ShapeFunction::NPOINTS
staticprotected

Definition at line 278 of file WellboreSimulatorFEM.h.

◆ jacobian_residual_size

template<typename ShapeFunction, int GlobalDim>
int const ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::jacobian_residual_size = 1
static

Definition at line 129 of file WellboreSimulatorFEM.h.

◆ pressure_index

template<typename ShapeFunction, int GlobalDim>
const int ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::pressure_index = 0
staticprotected

Definition at line 272 of file WellboreSimulatorFEM.h.

Referenced by assemble().

◆ pressure_size

template<typename ShapeFunction, int GlobalDim>
const int ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::pressure_size = ShapeFunction::NPOINTS
staticprotected

Definition at line 276 of file WellboreSimulatorFEM.h.

◆ velocity_index

template<typename ShapeFunction, int GlobalDim>
const int ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::velocity_index = ShapeFunction::NPOINTS
staticprotected

Definition at line 273 of file WellboreSimulatorFEM.h.

Referenced by assemble().

◆ velocity_size

template<typename ShapeFunction, int GlobalDim>
const int ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::velocity_size = ShapeFunction::NPOINTS
staticprotected

Definition at line 277 of file WellboreSimulatorFEM.h.


The documentation for this class was generated from the following files: