OGS
WellboreSimulatorFEM-impl.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <Eigen/Dense>
7#include <vector>
8
16
17namespace ProcessLib
18{
19namespace WellboreSimulator
20{
21template <typename ShapeFunction, int GlobalDim>
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)
26{
27 auto const local_matrix_size = local_x.size();
28
29 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
30
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);
37
38 // Get block matrices
39 auto Mvv = local_M.template block<velocity_size, velocity_size>(
41
42 auto Mhp = local_M.template block<enthalpy_size, pressure_size>(
44 auto Mhh = local_M.template block<enthalpy_size, enthalpy_size>(
46
47 auto Kpv = local_K.template block<pressure_size, velocity_size>(
49
50 auto Kvp = local_K.template block<velocity_size, pressure_size>(
52 auto Kvv = local_K.template block<velocity_size, velocity_size>(
54
55 auto Khh = local_K.template block<enthalpy_size, enthalpy_size>(
57
58 auto Bp = local_b.template segment<pressure_size>(pressure_index);
59 auto Bv = local_b.template segment<velocity_size>(velocity_index);
60 auto Bh = local_b.template segment<enthalpy_size>(enthalpy_index);
61
62 unsigned const n_integration_points =
63 _integration_method.getNumberOfPoints();
64
66 pos.setElementID(_element.getID());
67
68 auto const& b = _process_data.specific_body_force;
69
71
72 // Get material properties
73 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
74 auto const& liquid_phase =
76 auto const& gas_phase = medium.phase(MaterialPropertyLib::PhaseName::Gas);
77
78 // Get wellbore parameters
79 // casing thickness
80 auto const t_ca = _process_data.wellbore.casing_thickness(t, pos)[0];
81 // wellbore radius
82 auto const r_w = _process_data.wellbore.diameter(t, pos)[0] / 2;
83
84 // pipe thickness
85 auto const t_p = _process_data.wellbore.pipe_thickness(t, pos)[0];
86
87 // roughness of the wellbore
88 auto const xi = _process_data.wellbore.roughness(t, pos)[0];
89 // pipe outer radius
90 auto const r_o = r_w - t_ca;
91 // pipe inner radius
92 auto const r_i = r_o - t_p;
93
94 // get reservoir properties
95 NodalVectorType T_r =
96 _process_data.reservoir_properties.temperature.getNodalValuesOnElement(
97 _element, t);
98 NodalVectorType p_r =
99 _process_data.reservoir_properties.pressure.getNodalValuesOnElement(
100 _element, t);
102 _process_data.productivity_index.getNodalValuesOnElement(_element, t);
103 auto const k_r =
104 _process_data.reservoir_properties.thermal_conductivity(t, pos)[0];
105 auto const rho_r = _process_data.reservoir_properties.density(t, pos)[0];
106 auto const c_r =
107 _process_data.reservoir_properties.specific_heat_capacity(t, pos)[0];
108
109 for (unsigned ip(0); ip < n_integration_points; ip++)
110 {
111 auto& ip_data = _ip_data[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;
121
122 pos = {
123 std::nullopt, _element.getID(),
126 _element, N))};
127
128 double p_int_pt = 0.0;
129 double v_int_pt = 0.0;
130 double h_int_pt = 0.0;
131
132 NumLib::shapeFunctionInterpolate(local_x, N, p_int_pt, v_int_pt,
133 h_int_pt);
134
135 double p_prev_int_pt = 0.0;
136 double v_prev_int_pt = 0.0;
137 double h_prev_int_pt = 0.0;
138
139 NumLib::shapeFunctionInterpolate(local_x_prev, N, p_prev_int_pt,
140 v_prev_int_pt, h_prev_int_pt);
141
142 double vdot_int_pt = (v_int_pt - v_prev_int_pt) / dt;
143
144 // calculate fluid properties
145
146 const double pi = std::numbers::pi;
147
148 vars.liquid_phase_pressure = p_int_pt;
149 vars.enthalpy = h_int_pt;
150
151 double liquid_water_density =
152 liquid_phase
154 .template value<double>(vars, pos, t, dt);
155 double const vapour_water_density =
156 gas_phase
158 .template value<double>(vars, pos, t, dt);
159
160 double const h_sat_liq_w =
161 liquid_phase
162 .property(
164 .template value<double>(vars, pos, t, dt);
165 double const h_sat_vap_w =
166 gas_phase
167 .property(
169 .template value<double>(vars, pos, t, dt);
170
171 // TODO add a function to calculate dryness with constrain of 0
172 // to 1.
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;
176
177 double const T_int_pt =
178 (dryness == 0)
179 ? liquid_phase
181 .template value<double>(vars, pos, t, dt)
182 : gas_phase
184 saturation_temperature)
185 .template value<double>(vars, pos, t, dt);
186 temperature = T_int_pt;
187 vars.temperature = T_int_pt;
188
189 // For the calculation of the void fraction of vapour,
190 // see Rohuani, Z., and E. Axelsson. "Calculation of volume void
191 // fraction in a subcooled and quality region." International
192 // Journal of Heat and Mass Transfer 17 (1970): 383-393.
193
194 // profile parameter of drift flux
195 double C_0 = 1 + 0.12 * (1 - dryness);
196
197 // For the surface tension calculation, see
198 // Cooper, J. R., and R. B. Dooley. "IAPWS release on surface
199 // tension of ordinary water substance." International Association
200 // for the Properties of Water and Steam (1994).
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));
204 // drift flux velocity
205 double const u_gu =
206 1.18 * (1 - dryness) *
207 std::pow((9.81) * sigma_gl *
208 (liquid_water_density - vapour_water_density),
209 0.25) /
210 std::pow(liquid_water_density, 0.5);
211
212 // solving void fraction of vapor: Rouhani-Axelsson
213 double alpha = 0;
214 if (dryness != 0)
215 {
216 // Local Newton solver
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;
222
223 Eigen::PartialPivLU<LocalJacobianMatrix> linear_solver(1);
224
225 auto const update_residual = [&](LocalResidualVector& residual)
226 {
227 calculateResidual(alpha, vapour_water_density,
228 liquid_water_density, v_int_pt, dryness, C_0,
229 u_gu, residual);
230 };
231
232 auto const update_jacobian = [&](LocalJacobianMatrix& jacobian)
233 {
235 alpha, vapour_water_density, liquid_water_density, v_int_pt,
236 dryness, C_0, u_gu,
237 jacobian); // for solution dependent Jacobians
238 };
239
240 auto const update_solution =
241 [&](LocalUnknownVector const& increment)
242 {
243 // increment solution vectors
244 alpha += increment[0];
245 };
246
247 const int maximum_iterations(20);
248 const double residuum_tolerance(1.e-10);
249 const double increment_tolerance(0);
250
251 auto newton_solver = NumLib::NewtonRaphson(
252 linear_solver, update_jacobian, update_residual,
253 update_solution,
254 {maximum_iterations, residuum_tolerance, increment_tolerance});
255
256 auto const success_iterations = newton_solver.solve(J_loc);
257
258 if (!success_iterations)
259 {
260 WARN(
261 "Attention! Steam void fraction has not been correctly "
262 "calculated!");
263 }
264 }
265
266 vapor_volume_frac = alpha;
267
268 if (alpha == 0)
269 {
270 liquid_water_density =
271 liquid_phase
273 .template value<double>(vars, pos, t, dt);
274 }
275
276 mix_density =
277 vapour_water_density * alpha + liquid_water_density * (1 - alpha);
278
279 auto& mix_density_prev = ip_data.mix_density_prev;
280 vars.density = mix_density;
281
282 auto const rho_dot = (mix_density - mix_density_prev) / dt;
283
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 =
289 (alpha == 0) ? 0
290 : dryness * mix_density * v_int_pt /
291 (alpha * vapour_water_density);
292
293 vapor_mass_flowrate = vapor_water_velocity_act * vapour_water_density *
294 pi * r_i * r_i * alpha;
295
296 liquid_mass_flowrate = liquid_water_velocity_act *
297 liquid_water_density * pi * r_i * r_i *
298 (1 - alpha);
299
300 // Slip parameter between two phases,
301 // see Akbar, Somaieh, N. Fathianpour, and Rafid Al Khoury. "A finite
302 // element model for high enthalpy two-phase flow in geothermal
303 // wellbores." Renewable Energy 94 (2016): 223-236.
304 double const gamma =
305 alpha * liquid_water_density * vapour_water_density * mix_density /
306 (1 - alpha) /
307 std::pow((alpha * C_0 * vapour_water_density +
308 (1 - alpha * C_0) * liquid_water_density),
309 2) *
310 std::pow((C_0 - 1) * v_int_pt + u_gu, 2);
311
312 double const miu =
314 .template value<double>(vars, pos, t, dt);
315 double const Re = mix_density * v_int_pt * 2 * r_i / miu;
316
317 // Wall friction coefficient,
318 // Musivand Arzanfudi, Mehdi, and Rafid Al‐Khoury. "A compressible
319 // two‐fluid multiphase model for CO2 leakage through a wellbore."
320 // International Journal for Numerical Methods in Fluids 77.8 (2015):
321 // 477-507.
322 double f = 0.0;
323 if (Re > 10 && Re <= 2400)
324 {
325 f = 16 / Re;
326 }
327 else if (Re > 2400)
328 {
329 f = std::pow(std::log(xi / 3.7 / r_i) -
330 5.02 / Re * std::log(xi / 3.7 / r_i + 13 / Re),
331 -2) /
332 16;
333 }
334
335 double Q_hx = 0;
336 double const T_r_int_pt = N.dot(T_r);
337 // conductive heat exchange between wellbore and formation
338 if (_process_data.has_heat_exchange_with_formation)
339 {
340 // See Zhang, Pan, Pruess, Finsterle (2011). A time-convolution
341 // approach for modeling heat exchange between a wellbore and
342 // surrounding formation. Geothermics 40, 261-266.
343 const double alpha_r = k_r / rho_r / c_r;
344 const double t_d = alpha_r * t / (r_i * r_i);
345
346 double beta;
347 if (t_d < 2.8)
348 {
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;
351 }
352 else
353 {
354 beta = 2 * (1 / (std::log(4 * t_d) - 2 * 0.57722) -
355 0.57722 /
356 std::pow((std::log(4 * t_d) - 2 * 0.57722), 2));
357 }
358
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;
361 }
362
363 // mass exchange with reservoir
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);
367
368 // advective momentum and energy exchange with reservoir due to the mass
369 // exchange
370 double Q_mom = 0;
371 double Q_ene = 0;
372 if (Q_mx != 0)
373 {
374 Q_mom = Q_mx * v_int_pt;
375 // Only single-phase liquid condition is considered now
376 // TODO: update the two-phase flowing enthalpy from the feed zone.
377 vars.liquid_phase_pressure = p_r_int_pt;
378 vars.temperature = T_r_int_pt;
379 double const h_fres =
380 liquid_phase
382 .template value<double>(vars, pos, t, dt);
383 Q_ene = Q_mx * h_fres;
384 }
385
386 // M matrix assembly
387 Mvv.noalias() += w * N.transpose() * mix_density * N;
388
389 Mhp.noalias() += -w * N.transpose() * N;
390 Mhh.noalias() += w * N.transpose() * mix_density * N;
391
392 // K matrix assembly
393 Kpv.noalias() += w * dNdx.transpose() * N * mix_density;
394
395 Kvp.noalias() += w * N.transpose() * dNdx;
396 Kvv.noalias() += w * N.transpose() * rho_dot * N;
397
398 Khh.noalias() += w * N.transpose() * mix_density * v_int_pt * dNdx;
399
400 // b matrix assembly
401 Bp.noalias() += w * N.transpose() * rho_dot + w * N.transpose() * Q_mx;
402
403 Bv.noalias() +=
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;
409
410 Bh.noalias() +=
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 *
414 v_int_pt +
415 w * N.transpose() * (Q_hx / pi / r_i / r_i) -
416 w * N.transpose() * Q_ene;
417
418 if (_process_data.has_gravity)
419 {
420 NodalVectorType gravity_operator =
421 N.transpose() * b * w * _element_direction[2];
422
423 Bv.noalias() += gravity_operator * mix_density;
424 Bh.noalias() += gravity_operator * mix_density * v_int_pt;
425 }
426 }
427
428 // debugging
429 // std::string sep = "\n----------------------------------------\n";
430 // Eigen::IOFormat CleanFmt(6, 0, ", ", "\n", "[", "]");
431 // std::cout << local_M.format(CleanFmt) << sep;
432 // std::cout << local_K.format(CleanFmt) << sep;
433 // std::cout << local_b.format(CleanFmt) << sep;
434}
435
436} // namespace WellboreSimulator
437} // namespace ProcessLib
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:34
const double PI
Definition QArrow.h:9
std::optional< int > solve(JacobianMatrix &jacobian) const
void setElementID(std::size_t element_id)
typename ShapeMatricesType::NodalVectorType NodalVectorType
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)
NumLib::GenericIntegrationMethod const & _integration_method
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
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)
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
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 &)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)