www.mooseframework.org
PorousFlowPeacemanBorehole.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 #include "RotationMatrix.h"
12 #include "Function.h"
13 
15 
16 template <>
17 InputParameters
19 {
20  InputParameters params = validParams<PorousFlowLineSink>();
21  params.addRequiredParam<FunctionName>(
22  "character",
23  "If zero then borehole does nothing. If positive the borehole acts as a sink "
24  "(production well) for porepressure > borehole pressure, and does nothing "
25  "otherwise. If negative the borehole acts as a source (injection well) for "
26  "porepressure < borehole pressure, and does nothing otherwise. The flow rate "
27  "to/from the borehole is multiplied by |character|, so usually character = +/- "
28  "1, but you can specify other quantities to provide an overall scaling to the "
29  "flow if you like.");
30  params.addRequiredParam<Real>("bottom_p_or_t",
31  "For function_of=pressure, this parameter is the "
32  "pressure at the bottom of the borehole, "
33  "otherwise it is the temperature at the bottom of "
34  "the borehole");
35  params.addRequiredParam<RealVectorValue>(
36  "unit_weight",
37  "(fluid_density*gravitational_acceleration) as a vector pointing downwards. "
38  "Note that the borehole pressure at a given z position is bottom_p_or_t + "
39  "unit_weight*(q - q_bottom), where q=(x,y,z) and q_bottom=(x,y,z) of the "
40  "bottom point of the borehole. The analogous formula holds for "
41  "function_of=temperature. If you don't want bottomhole pressure (or "
42  "temperature) to vary in the borehole just set unit_weight=0. Typical value "
43  "is = (0,0,-1E4), for water");
44  params.addParam<Real>("re_constant",
45  0.28,
46  "The dimensionless constant used in evaluating the borehole effective "
47  "radius. This depends on the meshing scheme. Peacemann "
48  "finite-difference calculations give 0.28, while for rectangular finite "
49  "elements the result is closer to 0.1594. (See Eqn(4.13) of Z Chen, Y "
50  "Zhang, Well flow models for various numerical methods, Int J Num "
51  "Analysis and Modeling, 3 (2008) 375-388.)");
52  params.addParam<Real>("well_constant",
53  -1.0,
54  "Usually this is calculated internally from the element geometry, the "
55  "local borehole direction and segment length, and the permeability. "
56  "However, if this parameter is given as a positive number then this "
57  "number is used instead of the internal calculation. This speeds up "
58  "computation marginally. re_constant becomes irrelevant");
59  params.addClassDescription(
60  "Approximates a borehole in the mesh using the Peaceman approach, ie "
61  "using a number of point sinks with given radii whose positions are "
62  "read from a file. NOTE: if you are using PorousFlowPorosity that depends on volumetric "
63  "strain, you should set strain_at_nearest_qp=true in your GlobalParams, to ensure the nodal "
64  "Porosity Material uses the volumetric strain at the Dirac quadpoints, and can therefore be "
65  "computed");
66  return params;
67 }
68 
69 PorousFlowPeacemanBorehole::PorousFlowPeacemanBorehole(const InputParameters & parameters)
70  : PorousFlowLineSink(parameters),
71  _character(getFunction("character")),
72  _p_bot(getParam<Real>("bottom_p_or_t")),
73  _unit_weight(getParam<RealVectorValue>("unit_weight")),
74  _re_constant(getParam<Real>("re_constant")),
75  _well_constant(getParam<Real>("well_constant")),
76  _has_permeability(
77  hasMaterialProperty<RealTensorValue>("PorousFlow_permeability_qp") &&
78  hasMaterialProperty<std::vector<RealTensorValue>>("dPorousFlow_permeability_qp_dvar")),
79  _has_thermal_conductivity(
80  hasMaterialProperty<RealTensorValue>("PorousFlow_thermal_conductivity_qp") &&
81  hasMaterialProperty<std::vector<RealTensorValue>>(
82  "dPorousFlow_thermal_conductivity_qp_dvar")),
83  _perm_or_cond(_p_or_t == PorTchoice::pressure
84  ? getMaterialProperty<RealTensorValue>("PorousFlow_permeability_qp")
85  : getMaterialProperty<RealTensorValue>("PorousFlow_thermal_conductivity_qp")),
86  _dperm_or_cond_dvar(
87  _p_or_t == PorTchoice::pressure
88  ? getMaterialProperty<std::vector<RealTensorValue>>("dPorousFlow_permeability_qp_dvar")
89  : getMaterialProperty<std::vector<RealTensorValue>>(
90  "dPorousFlow_thermal_conductivity_qp_dvar"))
91 {
93  mooseError("PorousFlowPeacemanBorehole: You have specified function_of=porepressure, but you "
94  "do not have a quadpoint permeability material");
96  mooseError("PorousFlowPeacemanBorehole: You have specified function_of=temperature, but you do "
97  "not have a quadpoint thermal_conductivity material");
98 
99  // construct the rotation matrix needed to rotate the permeability
100  const unsigned int num_pts = _zs.size();
101  _rot_matrix.resize(std::max(num_pts - 1, (unsigned)1));
102  for (unsigned int i = 0; i + 1 < num_pts; ++i)
103  {
104  const RealVectorValue v2(_xs[i + 1] - _xs[i], _ys[i + 1] - _ys[i], _zs[i + 1] - _zs[i]);
105  _rot_matrix[i] = RotationMatrix::rotVecToZ(v2);
106  }
107  if (num_pts == (unsigned)1)
108  _rot_matrix[0] = RotationMatrix::rotVecToZ(_line_direction);
109 }
110 
111 Real
112 PorousFlowPeacemanBorehole::wellConstant(const RealTensorValue & perm,
113  const RealTensorValue & rot,
114  const Real & half_len,
115  const Elem * ele,
116  const Real & rad) const
117 // Peaceman's form for the borehole well constant
118 {
119  if (_well_constant > 0)
120  return _well_constant;
121 
122  // rot_perm has its "2" component lying along the half segment.
123  // We want to determine the eigenvectors of rot(0:1, 0:1), since, when
124  // rotated back to the original frame we will determine the element
125  // lengths along these directions
126  const RealTensorValue rot_perm = (rot * perm) * rot.transpose();
127  const Real trace2D = rot_perm(0, 0) + rot_perm(1, 1);
128  const Real det2D = rot_perm(0, 0) * rot_perm(1, 1) - rot_perm(0, 1) * rot_perm(1, 0);
129  const Real sq = std::sqrt(std::max(0.25 * trace2D * trace2D - det2D,
130  0.0)); // the std::max accounts for wierdo precision loss
131  const Real eig_val1 = 0.5 * trace2D + sq;
132  const Real eig_val2 = 0.5 * trace2D - sq;
133  RealVectorValue eig_vec1, eig_vec2;
134  if (sq > std::abs(trace2D) * 1E-7) // matrix is not a multiple of the identity (1E-7 accounts for
135  // precision in a crude way)
136  {
137  if (rot_perm(1, 0) != 0)
138  {
139  eig_vec1(0) = eig_val1 - rot_perm(1, 1);
140  eig_vec1(1) = rot_perm(1, 0);
141  eig_vec2(0) = eig_val2 - rot_perm(1, 1);
142  eig_vec2(1) = rot_perm(1, 0);
143  }
144  else if (rot_perm(0, 1) != 0)
145  {
146  eig_vec1(0) = rot_perm(0, 1);
147  eig_vec1(1) = eig_val1 - rot_perm(0, 0);
148  eig_vec2(0) = rot_perm(0, 1);
149  eig_vec2(1) = eig_val2 - rot_perm(0, 0);
150  }
151  else // off diagonal terms are both zero
152  {
153  eig_vec1(0) = 1.0;
154  eig_vec2(1) = 1.0;
155  }
156  }
157  else // matrix is basically a multiple of the identity
158  {
159  eig_vec1(0) = 1.0;
160  eig_vec2(1) = 1.0;
161  }
162 
163  // finally, rotate these to original frame and normalise
164  eig_vec1 = rot.transpose() * eig_vec1;
165  eig_vec1 /= std::sqrt(eig_vec1 * eig_vec1);
166  eig_vec2 = rot.transpose() * eig_vec2;
167  eig_vec2 /= std::sqrt(eig_vec2 * eig_vec2);
168 
169  // find the "length" of the element in these directions
170  // TODO - maybe better to use variance than max&min
171  Real max1 = eig_vec1 * ele->point(0);
172  Real max2 = eig_vec2 * ele->point(0);
173  Real min1 = max1;
174  Real min2 = max2;
175  Real proj;
176  for (unsigned int i = 1; i < ele->n_nodes(); i++)
177  {
178  proj = eig_vec1 * ele->point(i);
179  max1 = (max1 < proj) ? proj : max1;
180  min1 = (min1 < proj) ? min1 : proj;
181 
182  proj = eig_vec2 * ele->point(i);
183  max2 = (max2 < proj) ? proj : max2;
184  min2 = (min2 < proj) ? min2 : proj;
185  }
186  const Real ll1 = max1 - min1;
187  const Real ll2 = max2 - min2;
188 
189  Real r0;
190  if (eig_val1 <= 0.0)
191  r0 = _re_constant * ll1;
192  else if (eig_val2 <= 0.0)
193  r0 = _re_constant * ll2;
194  else
195  r0 = _re_constant * std::sqrt(std::sqrt(eig_val1 / eig_val2) * std::pow(ll2, 2) +
196  std::sqrt(eig_val2 / eig_val1) * std::pow(ll1, 2)) /
197  (std::pow(eig_val1 / eig_val2, 0.25) + std::pow(eig_val2 / eig_val1, 0.25));
198 
199  const Real effective_perm = (det2D >= 0.0 ? std::sqrt(det2D) : 0.0);
200 
201  const Real halfPi = acos(0.0);
202 
203  if (r0 <= rad)
204  mooseError("The effective element size (about 0.2-times-true-ele-size) for an element "
205  "containing a Peaceman-type borehole must be (much) larger than the borehole radius "
206  "for the Peaceman formulation to be correct. Your element has effective size ",
207  r0,
208  " and the borehole radius is ",
209  rad,
210  "\n");
211 
212  return 4 * halfPi * effective_perm * half_len / std::log(r0 / rad);
213 }
214 
215 Real
216 PorousFlowPeacemanBorehole::computeQpBaseOutflow(unsigned current_dirac_ptid) const
217 {
218  const Real character = _character.value(_t, _q_point[_qp]);
219  if (character == 0.0)
220  return 0.0;
221 
222  const Real bh_pressure = _p_bot + _unit_weight * (_q_point[_qp] - _bottom_point);
223  const Real pp = ptqp();
224 
225  Real outflow = 0.0; // this is the flow rate from porespace out of the system
226 
227  if (current_dirac_ptid > 0)
228  // contribution from half-segment "behind" this point (must have >1 point for
229  // current_dirac_ptid>0)
230  {
231  if ((character < 0.0 && pp < bh_pressure) || (character > 0.0 && pp > bh_pressure))
232  {
233  // injection, so outflow<0 || production, so outflow>0
234  const Real wc = wellConstant(_perm_or_cond[_qp],
235  _rot_matrix[current_dirac_ptid - 1],
236  _half_seg_len[current_dirac_ptid - 1],
237  _current_elem,
238  _rs[current_dirac_ptid]);
239  outflow += wc * (pp - bh_pressure);
240  }
241  }
242 
243  if (current_dirac_ptid + 1 < _zs.size() || _zs.size() == 1)
244  // contribution from half-segment "ahead of" this point, or we only have one point
245  {
246  if ((character < 0.0 && pp < bh_pressure) || (character > 0.0 && pp > bh_pressure))
247  {
248  // injection, so outflow<0 || // production, so outflow>0
249  const Real wc = wellConstant(_perm_or_cond[_qp],
250  _rot_matrix[current_dirac_ptid],
251  _half_seg_len[current_dirac_ptid],
252  _current_elem,
253  _rs[current_dirac_ptid]);
254  outflow += wc * (pp - bh_pressure);
255  }
256  }
257 
258  return outflow * _test[_i][_qp] * std::abs(character);
259 }
260 
261 void
263  unsigned current_dirac_ptid,
264  Real & outflow,
265  Real & outflowp) const
266 {
267  outflow = 0.0;
268  outflowp = 0.0;
269 
270  const Real character = _character.value(_t, _q_point[_qp]);
271  if (character == 0.0)
272  return;
273 
275  return;
276  const unsigned pvar = _dictator.porousFlowVariableNum(jvar);
277 
278  const Real bh_pressure = _p_bot + _unit_weight * (_q_point[_qp] - _bottom_point);
279  const Real pp = ptqp();
280  const Real pp_prime = dptqp(pvar) * _phi[_j][_qp];
281 
282  if (current_dirac_ptid > 0)
283  // contribution from half-segment "behind" this point
284  {
285  if ((character < 0.0 && pp < bh_pressure) || (character > 0.0 && pp > bh_pressure))
286  {
287  // injection, so outflow<0 || // production, so outflow>0
288  const Real wc = wellConstant(_perm_or_cond[_qp],
289  _rot_matrix[current_dirac_ptid - 1],
290  _half_seg_len[current_dirac_ptid - 1],
291  _current_elem,
292  _rs[current_dirac_ptid]);
293  outflowp += wc * pp_prime;
294  outflow += wc * (pp - bh_pressure);
295  }
296  }
297 
298  if (current_dirac_ptid < _zs.size() - 1 || _zs.size() == 1)
299  // contribution from half-segment "ahead of" this point
300  {
301  if ((character < 0.0 && pp < bh_pressure) || (character > 0.0 && pp > bh_pressure))
302  {
303  // injection, so outflow<0 || // production, so outflow>0
304  const Real wc = wellConstant(_perm_or_cond[_qp],
305  _rot_matrix[current_dirac_ptid],
306  _half_seg_len[current_dirac_ptid],
307  _current_elem,
308  _rs[current_dirac_ptid]);
309  outflowp += wc * pp_prime;
310  outflow += wc * (pp - bh_pressure);
311  }
312  }
313 
314  outflowp *= _test[_i][_qp] * std::abs(character);
315  outflow *= _test[_i][_qp] * std::abs(character);
316 }
PorousFlowLineSink
Approximates a line sink a sequence of Dirac Points.
Definition: PorousFlowLineSink.h:24
PorousFlowPeacemanBorehole::_has_thermal_conductivity
const bool _has_thermal_conductivity
Whether there is a quadpoint thermal conductivity material (for error checking)
Definition: PorousFlowPeacemanBorehole.h:61
pow
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
Definition: ExpressionBuilder.h:673
PorousFlowPeacemanBorehole::PorousFlowPeacemanBorehole
PorousFlowPeacemanBorehole(const InputParameters &parameters)
Creates a new PorousFlowPeacemanBorehole This reads the file containing the lines of the form radius ...
Definition: PorousFlowPeacemanBorehole.C:69
PorousFlowLineSink::PorTchoice
PorTchoice
whether the flux is a function of pressure or temperature
Definition: PorousFlowLineSink.h:91
PorousFlowPeacemanBorehole
Approximates a borehole by a sequence of Dirac Points.
Definition: PorousFlowPeacemanBorehole.h:22
PorousFlowPeacemanBorehole::_has_permeability
const bool _has_permeability
Whether there is a quadpoint permeability material (for error checking)
Definition: PorousFlowPeacemanBorehole.h:58
validParams< PorousFlowLineSink >
InputParameters validParams< PorousFlowLineSink >()
Definition: PorousFlowLineSink.C:15
PorousFlowLineGeometry::_xs
std::vector< Real > _xs
x points of the borehole
Definition: PorousFlowLineGeometry.h:52
PorousFlowLineGeometry::_rs
std::vector< Real > _rs
Radii of the borehole.
Definition: PorousFlowLineGeometry.h:49
PorousFlowDictator::notPorousFlowVariable
bool notPorousFlowVariable(unsigned int moose_var_num) const
Returns true if moose_var_num is not a porous flow variabe.
Definition: PorousFlowDictator.C:161
PorousFlowPeacemanBorehole::computeQpBaseOutflowJacobian
void computeQpBaseOutflowJacobian(unsigned jvar, unsigned current_dirac_ptid, Real &outflow, Real &outflowp) const override
Calculates the BaseOutflow as well as its derivative wrt jvar. Derived classes should override this.
Definition: PorousFlowPeacemanBorehole.C:262
PorousFlowLineSink::PorTchoice::pressure
PorousFlowLineSink::PorTchoice::temperature
PorousFlowPeacemanBorehole::_perm_or_cond
const MaterialProperty< RealTensorValue > & _perm_or_cond
Permeability or conductivity of porous material.
Definition: PorousFlowPeacemanBorehole.h:64
PorousFlowPeacemanBorehole.h
PorousFlowDictator::porousFlowVariableNum
unsigned int porousFlowVariableNum(unsigned int moose_var_num) const
The PorousFlow variable number.
Definition: PorousFlowDictator.C:135
PorousFlowPeacemanBorehole::_character
const Function & _character
If positive then the borehole acts as a sink (producion well) for porepressure > borehole pressure,...
Definition: PorousFlowPeacemanBorehole.h:43
PorousFlowPeacemanBorehole::wellConstant
Real wellConstant(const RealTensorValue &perm, const RealTensorValue &rot, const Real &half_len, const Elem *ele, const Real &rad) const
Calculates Peaceman's form of the borehole well constant Z Chen, Y Zhang, Well flow models for variou...
Definition: PorousFlowPeacemanBorehole.C:112
PorousFlowPeacemanBorehole::_unit_weight
const RealVectorValue _unit_weight
Unit weight of fluid in borehole (for calculating bottomhole pressure at each Dirac Point)
Definition: PorousFlowPeacemanBorehole.h:49
PorousFlowLineGeometry::_bottom_point
Point _bottom_point
The bottom point of the borehole (where bottom_pressure is defined)
Definition: PorousFlowLineGeometry.h:61
PorousFlowPeacemanBorehole::computeQpBaseOutflow
Real computeQpBaseOutflow(unsigned current_dirac_ptid) const override
Returns the flux from the line sink (before modification by mobility, etc). Derived classes should ov...
Definition: PorousFlowPeacemanBorehole.C:216
PorousFlowPeacemanBorehole::_p_bot
const Real _p_bot
Bottomhole pressure of borehole.
Definition: PorousFlowPeacemanBorehole.h:46
validParams< PorousFlowPeacemanBorehole >
InputParameters validParams< PorousFlowPeacemanBorehole >()
Definition: PorousFlowPeacemanBorehole.C:18
registerMooseObject
registerMooseObject("PorousFlowApp", PorousFlowPeacemanBorehole)
PorousFlowLineSink::_p_or_t
enum PorousFlowLineSink::PorTchoice _p_or_t
PorousFlowLineSink::dptqp
Real dptqp(unsigned pvar) const
If _p_or_t==0, then returns d(quadpoint porepressure)/d(PorousFlow variable), else returns d(quadpoin...
Definition: PorousFlowLineSink.C:341
PorousFlowLineGeometry::_zs
std::vector< Real > _zs
z points of borehole
Definition: PorousFlowLineGeometry.h:58
PorousFlowLineGeometry::_half_seg_len
std::vector< Real > _half_seg_len
0.5*(length of polyline segments between points)
Definition: PorousFlowLineGeometry.h:64
PorousFlowLineGeometry::_line_direction
const RealVectorValue _line_direction
Line direction. This is only used if there is only one borehole point.
Definition: PorousFlowLineGeometry.h:39
PorousFlowLineGeometry::_ys
std::vector< Real > _ys
y points of the borehole
Definition: PorousFlowLineGeometry.h:55
PorousFlowLineSink::ptqp
Real ptqp() const
If _p_or_t==0, then returns the quadpoint porepressure, else returns the quadpoint temperature.
Definition: PorousFlowLineSink.C:335
PorousFlowLineSink::_dictator
const PorousFlowDictator & _dictator
PorousFlowDictator UserObject.
Definition: PorousFlowLineSink.h:60
NS::pressure
const std::string pressure
Definition: NS.h:25
PorousFlowPeacemanBorehole::_rot_matrix
std::vector< RealTensorValue > _rot_matrix
Rotation matrix used in well_constant calculation.
Definition: PorousFlowPeacemanBorehole.h:70