23 double const t,
double const dt, std::vector<double>
const& local_x,
24 std::vector<double>
const& local_x_prev, std::vector<double>& local_M_data,
25 std::vector<double>& local_K_data, std::vector<double>& local_b_data)
27 auto const local_matrix_size = local_x.size();
29 assert(local_matrix_size == ShapeFunction::NPOINTS *
NUM_NODAL_DOF);
32 local_M_data, local_matrix_size, local_matrix_size);
34 local_K_data, local_matrix_size, local_matrix_size);
36 local_b_data, local_matrix_size);
39 auto Mvv = local_M.template block<velocity_size, velocity_size>(
42 auto Mhp = local_M.template block<enthalpy_size, pressure_size>(
44 auto Mhh = local_M.template block<enthalpy_size, enthalpy_size>(
47 auto Kpv = local_K.template block<pressure_size, velocity_size>(
50 auto Kvp = local_K.template block<velocity_size, pressure_size>(
52 auto Kvv = local_K.template block<velocity_size, velocity_size>(
55 auto Khh = local_K.template block<enthalpy_size, enthalpy_size>(
62 unsigned const n_integration_points =
74 auto const& liquid_phase =
80 auto const t_ca =
_process_data.wellbore.casing_thickness(t, pos)[0];
82 auto const r_w =
_process_data.wellbore.diameter(t, pos)[0] / 2;
85 auto const t_p =
_process_data.wellbore.pipe_thickness(t, pos)[0];
90 auto const r_o = r_w - t_ca;
92 auto const r_i = r_o - t_p;
96 _process_data.reservoir_properties.temperature.getNodalValuesOnElement(
99 _process_data.reservoir_properties.pressure.getNodalValuesOnElement(
104 _process_data.reservoir_properties.thermal_conductivity(t, pos)[0];
105 auto const rho_r =
_process_data.reservoir_properties.density(t, pos)[0];
107 _process_data.reservoir_properties.specific_heat_capacity(t, pos)[0];
109 for (
unsigned ip(0); ip < n_integration_points; ip++)
112 auto const& N = ip_data.N;
113 auto const& dNdx = ip_data.dNdx;
114 auto const& w = ip_data.integration_weight;
115 auto& mix_density = ip_data.mix_density;
116 auto& temperature = ip_data.temperature;
117 auto& steam_mass_frac = ip_data.dryness;
118 auto& vapor_volume_frac = ip_data.vapor_volume_fraction;
119 auto& vapor_mass_flowrate = ip_data.vapor_mass_flow_rate;
120 auto& liquid_mass_flowrate = ip_data.liquid_mass_flow_rate;
128 double p_int_pt = 0.0;
129 double v_int_pt = 0.0;
130 double h_int_pt = 0.0;
135 double p_prev_int_pt = 0.0;
136 double v_prev_int_pt = 0.0;
137 double h_prev_int_pt = 0.0;
140 v_prev_int_pt, h_prev_int_pt);
142 double vdot_int_pt = (v_int_pt - v_prev_int_pt) / dt;
146 const double pi = std::numbers::pi;
151 double liquid_water_density =
154 .template value<double>(vars, pos, t, dt);
155 double const vapour_water_density =
158 .template value<double>(vars, pos, t, dt);
160 double const h_sat_liq_w =
164 .template value<double>(vars, pos, t, dt);
165 double const h_sat_vap_w =
169 .template value<double>(vars, pos, t, dt);
173 double const dryness = std::max(
174 0., (h_int_pt - h_sat_liq_w) / (h_sat_vap_w - h_sat_liq_w));
175 steam_mass_frac = dryness;
177 double const T_int_pt =
181 .template value<double>(vars, pos, t, dt)
184 saturation_temperature)
185 .template value<double>(vars, pos, t, dt);
186 temperature = T_int_pt;
195 double C_0 = 1 + 0.12 * (1 - dryness);
201 double const sigma_gl = 0.2358 *
202 std::pow((1 - T_int_pt / 647.096), 1.256) *
203 (1 - 0.625 * (1 - T_int_pt / 647.096));
206 1.18 * (1 - dryness) *
207 std::pow((9.81) * sigma_gl *
208 (liquid_water_density - vapour_water_density),
210 std::pow(liquid_water_density, 0.5);
217 using LocalJacobianMatrix =
218 Eigen::Matrix<double, 1, 1, Eigen::RowMajor>;
219 using LocalResidualVector = Eigen::Matrix<double, 1, 1>;
220 using LocalUnknownVector = Eigen::Matrix<double, 1, 1>;
221 LocalJacobianMatrix J_loc;
223 Eigen::PartialPivLU<LocalJacobianMatrix> linear_solver(1);
225 auto const update_residual = [&](LocalResidualVector& residual)
228 liquid_water_density, v_int_pt, dryness, C_0,
232 auto const update_jacobian = [&](LocalJacobianMatrix& jacobian)
235 alpha, vapour_water_density, liquid_water_density, v_int_pt,
240 auto const update_solution =
241 [&](LocalUnknownVector
const& increment)
244 alpha += increment[0];
247 const int maximum_iterations(20);
248 const double residuum_tolerance(1.e-10);
249 const double increment_tolerance(0);
252 linear_solver, update_jacobian, update_residual,
254 {maximum_iterations, residuum_tolerance, increment_tolerance});
256 auto const success_iterations = newton_solver.
solve(J_loc);
258 if (!success_iterations)
261 "Attention! Steam void fraction has not been correctly "
266 vapor_volume_frac = alpha;
270 liquid_water_density =
273 .template value<double>(vars, pos, t, dt);
277 vapour_water_density * alpha + liquid_water_density * (1 - alpha);
279 auto& mix_density_prev = ip_data.mix_density_prev;
282 auto const rho_dot = (mix_density - mix_density_prev) / dt;
284 double const liquid_water_velocity_act =
285 (alpha == 0) ? v_int_pt
286 : (1 - dryness) * mix_density * v_int_pt /
287 (1 - alpha) / liquid_water_density;
288 double const vapor_water_velocity_act =
290 : dryness * mix_density * v_int_pt /
291 (alpha * vapour_water_density);
293 vapor_mass_flowrate = vapor_water_velocity_act * vapour_water_density *
294 pi * r_i * r_i * alpha;
296 liquid_mass_flowrate = liquid_water_velocity_act *
297 liquid_water_density * pi * r_i * r_i *
305 alpha * liquid_water_density * vapour_water_density * mix_density /
307 std::pow((alpha * C_0 * vapour_water_density +
308 (1 - alpha * C_0) * liquid_water_density),
310 std::pow((C_0 - 1) * v_int_pt + u_gu, 2);
314 .template value<double>(vars, pos, t, dt);
315 double const Re = mix_density * v_int_pt * 2 * r_i / miu;
323 if (Re > 10 && Re <= 2400)
329 f = std::pow(std::log(xi / 3.7 / r_i) -
330 5.02 / Re * std::log(xi / 3.7 / r_i + 13 / Re),
336 double const T_r_int_pt = N.dot(T_r);
343 const double alpha_r = k_r / rho_r / c_r;
344 const double t_d = alpha_r * t / (r_i * r_i);
349 beta = std::pow((pi * t_d), -0.5) + 0.5 -
350 0.25 * std::pow((t_d / pi), 0.5) + 0.125 * t_d;
354 beta = 2 * (1 / (std::log(4 * t_d) - 2 * 0.57722) -
356 std::pow((std::log(4 * t_d) - 2 * 0.57722), 2));
359 const double P_c = 2 * pi * r_i;
360 Q_hx = P_c * k_r * (T_r_int_pt - T_int_pt) / r_i * beta;
364 double const p_r_int_pt = N.dot(p_r);
365 double const PI_int_pt = N.dot(
PI);
366 double Q_mx = PI_int_pt * (p_int_pt - p_r_int_pt);
374 Q_mom = Q_mx * v_int_pt;
379 double const h_fres =
382 .template value<double>(vars, pos, t, dt);
383 Q_ene = Q_mx * h_fres;
387 Mvv.noalias() += w * N.transpose() * mix_density * N;
389 Mhp.noalias() += -w * N.transpose() * N;
390 Mhh.noalias() += w * N.transpose() * mix_density * N;
393 Kpv.noalias() += w * dNdx.transpose() * N * mix_density;
395 Kvp.noalias() += w * N.transpose() * dNdx;
396 Kvv.noalias() += w * N.transpose() * rho_dot * N;
398 Khh.noalias() += w * N.transpose() * mix_density * v_int_pt * dNdx;
401 Bp.noalias() += w * N.transpose() * rho_dot + w * N.transpose() * Q_mx;
404 w * dNdx.transpose() * mix_density * v_int_pt * v_int_pt +
405 w * dNdx.transpose() * gamma -
406 w * N.transpose() * f * mix_density * std::abs(v_int_pt) *
407 v_int_pt / (4 * r_i) -
408 w * N.transpose() * Q_mom;
411 -1 / 2 * w * N.transpose() * rho_dot * v_int_pt * v_int_pt -
412 w * N.transpose() * mix_density * v_int_pt * vdot_int_pt +
413 1 / 2 * w * dNdx.transpose() * mix_density * v_int_pt * v_int_pt *
415 w * N.transpose() * (Q_hx / pi / r_i / r_i) -
416 w * N.transpose() * Q_ene;
423 Bv.noalias() += gravity_operator * mix_density;
424 Bh.noalias() += gravity_operator * mix_density * v_int_pt;