www.mooseframework.org
PeacemanBorehole.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 
10 #include "PeacemanBorehole.h"
11 #include "RotationMatrix.h"
12 
13 #include <fstream>
14 
15 template <>
16 InputParameters
18 {
19  InputParameters params = validParams<DiracKernel>();
20  params.addRequiredParam<FunctionName>(
21  "character",
22  "If zero then borehole does nothing. If positive the borehole acts as a sink "
23  "(production well) for porepressure > borehole pressure, and does nothing "
24  "otherwise. If negative the borehole acts as a source (injection well) for "
25  "porepressure < borehole pressure, and does nothing otherwise. The flow rate "
26  "to/from the borehole is multiplied by |character|, so usually character = +/- "
27  "1, but you can specify other quantities to provide an overall scaling to the "
28  "flow if you like.");
29  params.addRequiredParam<Real>("bottom_pressure", "Pressure at the bottom of the borehole");
30  params.addRequiredParam<RealVectorValue>(
31  "unit_weight",
32  "(fluid_density*gravitational_acceleration) as a vector pointing downwards. "
33  "Note that the borehole pressure at a given z position is bottom_pressure + "
34  "unit_weight*(p - p_bottom), where p=(x,y,z) and p_bottom=(x,y,z) of the "
35  "bottom point of the borehole. If you don't want bottomhole pressure to vary "
36  "in the borehole just set unit_weight=0. Typical value is = (0,0,-1E4)");
37  params.addRequiredParam<std::string>(
38  "point_file",
39  "The file containing the borehole radii and coordinates of the point sinks "
40  "that approximate the borehole. Each line in the file must contain a "
41  "space-separated radius and coordinate. Ie r x y z. The last point in the "
42  "file is defined as the borehole bottom, where the borehole pressure is "
43  "bottom_pressure. If your file contains just one point, you must also specify "
44  "the borehole_length and borehole_direction. Note that you will get "
45  "segementation faults if your points do not lie within your mesh!");
46  params.addRequiredParam<UserObjectName>(
47  "SumQuantityUO",
48  "User Object of type=RichardsSumQuantity in which to place the total "
49  "outflow from the borehole for each time step.");
50  params.addParam<Real>("re_constant",
51  0.28,
52  "The dimensionless constant used in evaluating the borehole effective "
53  "radius. This depends on the meshing scheme. Peacemann "
54  "finite-difference calculations give 0.28, while for rectangular finite "
55  "elements the result is closer to 0.1594. (See Eqn(4.13) of Z Chen, Y "
56  "Zhang, Well flow models for various numerical methods, Int J Num "
57  "Analysis and Modeling, 3 (2008) 375-388.)");
58  params.addParam<Real>("well_constant",
59  -1.0,
60  "Usually this is calculated internally from the element geometry, the "
61  "local borehole direction and segment length, and the permeability. "
62  "However, if this parameter is given as a positive number then this "
63  "number is used instead of the internal calculation. This speeds up "
64  "computation marginally. re_constant becomes irrelevant");
65  params.addRangeCheckedParam<Real>(
66  "borehole_length",
67  0.0,
68  "borehole_length>=0",
69  "Borehole length. Note this is only used if there is only one point in the point_file.");
70  params.addParam<RealVectorValue>(
71  "borehole_direction",
72  RealVectorValue(0, 0, 1),
73  "Borehole direction. Note this is only used if there is only one point in the point_file.");
74  params.addClassDescription("Approximates a borehole in the mesh using the Peaceman approach, ie "
75  "using a number of point sinks with given radii whose positions are "
76  "read from a file");
77  return params;
78 }
79 
80 PeacemanBorehole::PeacemanBorehole(const InputParameters & parameters)
81  : DiracKernel(parameters),
82  _re_constant(getParam<Real>("re_constant")),
83  _well_constant(getParam<Real>("well_constant")),
84  _borehole_length(getParam<Real>("borehole_length")),
85  _borehole_direction(getParam<RealVectorValue>("borehole_direction")),
86  _point_file(getParam<std::string>("point_file")),
87  _character(getFunction("character")),
88  _p_bot(getParam<Real>("bottom_pressure")),
89  _unit_weight(getParam<RealVectorValue>("unit_weight")),
90  _total_outflow_mass(
91  const_cast<RichardsSumQuantity &>(getUserObject<RichardsSumQuantity>("SumQuantityUO")))
92 {
93  // zero the outflow mass
95 
96  // open file
97  std::ifstream file(_point_file.c_str());
98  if (!file.good())
99  mooseError("Error opening file '" + _point_file + "' from a Peaceman-type Borehole.");
100 
101  // construct the arrays of radius, x, y and z
102  std::vector<Real> scratch;
103  while (parseNextLineReals(file, scratch))
104  {
105  if (scratch.size() >= 2)
106  {
107  _rs.push_back(scratch[0]);
108  _xs.push_back(scratch[1]);
109  if (scratch.size() >= 3)
110  _ys.push_back(scratch[2]);
111  else
112  _ys.push_back(0.0);
113  if (scratch.size() >= 4)
114  _zs.push_back(scratch[3]);
115  else
116  _zs.push_back(0.0);
117  }
118  }
119 
120  file.close();
121 
122  const int num_pts = _zs.size();
123  _bottom_point(0) = _xs[num_pts - 1];
124  _bottom_point(1) = _ys[num_pts - 1];
125  _bottom_point(2) = _zs[num_pts - 1];
126 
127  // construct the line-segment lengths between each point
128  _half_seg_len.resize(std::max(num_pts - 1, 1));
129  for (unsigned int i = 0; i + 1 < _xs.size(); ++i)
130  {
131  _half_seg_len[i] =
132  0.5 * std::sqrt(std::pow(_xs[i + 1] - _xs[i], 2) + std::pow(_ys[i + 1] - _ys[i], 2) +
133  std::pow(_zs[i + 1] - _zs[i], 2));
134  if (_half_seg_len[i] == 0)
135  mooseError("Peaceman-type borehole has a zero-segment length at (x,y,z) = ",
136  _xs[i],
137  " ",
138  _ys[i],
139  " ",
140  _zs[i],
141  "\n");
142  }
143  if (num_pts == 1)
145 
146  // construct the rotation matrix needed to rotate the permeability
147  _rot_matrix.resize(std::max(num_pts - 1, 1));
148  for (unsigned int i = 0; i + 1 < _xs.size(); ++i)
149  {
150  const RealVectorValue v2(_xs[i + 1] - _xs[i], _ys[i + 1] - _ys[i], _zs[i + 1] - _zs[i]);
151  _rot_matrix[i] = RotationMatrix::rotVecToZ(v2);
152  }
153  if (num_pts == 1)
154  _rot_matrix[0] = RotationMatrix::rotVecToZ(_borehole_direction);
155 }
156 
157 bool
158 PeacemanBorehole::parseNextLineReals(std::ifstream & ifs, std::vector<Real> & myvec)
159 // reads a space-separated line of floats from ifs and puts in myvec
160 {
161  std::string line;
162  myvec.clear();
163  bool gotline(false);
164  if (getline(ifs, line))
165  {
166  gotline = true;
167 
168  // Harvest floats separated by whitespace
169  std::istringstream iss(line);
170  Real f;
171  while (iss >> f)
172  {
173  myvec.push_back(f);
174  }
175  }
176  return gotline;
177 }
178 
179 void
181 {
182  // This function gets called just before the DiracKernel is evaluated
183  // so this is a handy place to zero this out.
185 
186  // Add point using the unique ID "i", let the DiracKernel take
187  // care of the caching. This should be fast after the first call,
188  // as long as the points don't move around.
189  for (unsigned int i = 0; i < _zs.size(); i++)
190  addPoint(Point(_xs[i], _ys[i], _zs[i]), i);
191 }
192 
193 Real
194 PeacemanBorehole::wellConstant(const RealTensorValue & perm,
195  const RealTensorValue & rot,
196  const Real & half_len,
197  const Elem * ele,
198  const Real & rad)
199 // Peaceman's form for the borehole well constant
200 {
201  if (_well_constant > 0)
202  return _well_constant;
203 
204  // rot_perm has its "2" component lying along the half segment
205  // we want to determine the eigenvectors of rot(0:1, 0:1), since, when
206  // rotated back to the original frame we will determine the element
207  // lengths along these directions
208  const RealTensorValue rot_perm = (rot * perm) * rot.transpose();
209  const Real trace2D = rot_perm(0, 0) + rot_perm(1, 1);
210  const Real det2D = rot_perm(0, 0) * rot_perm(1, 1) - rot_perm(0, 1) * rot_perm(1, 0);
211  const Real sq = std::sqrt(std::max(0.25 * trace2D * trace2D - det2D,
212  0.0)); // the std::max accounts for wierdo precision loss
213  const Real eig_val1 = 0.5 * trace2D + sq;
214  const Real eig_val2 = 0.5 * trace2D - sq;
215  RealVectorValue eig_vec1, eig_vec2;
216  if (sq > std::abs(trace2D) * 1E-7) // matrix is not a multiple of the identity (1E-7 accounts for
217  // precision in a crude way)
218  {
219  if (rot_perm(1, 0) != 0)
220  {
221  eig_vec1(0) = eig_val1 - rot_perm(1, 1);
222  eig_vec1(1) = rot_perm(1, 0);
223  eig_vec2(0) = eig_val2 - rot_perm(1, 1);
224  eig_vec2(1) = rot_perm(1, 0);
225  }
226  else if (rot_perm(0, 1) != 0)
227  {
228  eig_vec1(0) = rot_perm(0, 1);
229  eig_vec1(1) = eig_val1 - rot_perm(0, 0);
230  eig_vec2(0) = rot_perm(0, 1);
231  eig_vec2(1) = eig_val2 - rot_perm(0, 0);
232  }
233  else // off diagonal terms are both zero
234  {
235  eig_vec1(0) = 1;
236  eig_vec2(1) = 1;
237  }
238  }
239  else // matrix is basically a multiple of the identity
240  {
241  eig_vec1(0) = 1;
242  eig_vec2(1) = 1;
243  }
244 
245  // finally, rotate these to original frame and normalise
246  eig_vec1 = rot.transpose() * eig_vec1;
247  eig_vec1 /= std::sqrt(eig_vec1 * eig_vec1);
248  eig_vec2 = rot.transpose() * eig_vec2;
249  eig_vec2 /= std::sqrt(eig_vec2 * eig_vec2);
250 
251  // find the "length" of the element in these directions
252  // TODO - probably better to use variance than max&min
253  Real max1 = eig_vec1 * ele->point(0);
254  Real max2 = eig_vec2 * ele->point(0);
255  Real min1 = max1;
256  Real min2 = max2;
257  Real proj;
258  for (unsigned int i = 1; i < ele->n_nodes(); i++)
259  {
260  proj = eig_vec1 * ele->point(i);
261  max1 = (max1 < proj) ? proj : max1;
262  min1 = (min1 < proj) ? min1 : proj;
263 
264  proj = eig_vec2 * ele->point(i);
265  max2 = (max2 < proj) ? proj : max2;
266  min2 = (min2 < proj) ? min2 : proj;
267  }
268  const Real ll1 = max1 - min1;
269  const Real ll2 = max2 - min2;
270 
271  const Real r0 = _re_constant * std::sqrt(std::sqrt(eig_val1 / eig_val2) * std::pow(ll2, 2) +
272  std::sqrt(eig_val2 / eig_val1) * std::pow(ll1, 2)) /
273  (std::pow(eig_val1 / eig_val2, 0.25) + std::pow(eig_val2 / eig_val1, 0.25));
274 
275  const Real effective_perm = std::sqrt(det2D);
276 
277  const Real halfPi = acos(0.0);
278 
279  if (r0 <= rad)
280  mooseError("The effective element size (about 0.2-times-true-ele-size) for an element "
281  "containing a Peaceman-type borehole must be (much) larger than the borehole radius "
282  "for the Peaceman formulation to be correct. Your element has effective size ",
283  r0,
284  " and the borehole radius is ",
285  rad,
286  "\n");
287 
288  return 4 * halfPi * effective_perm * half_len / std::log(r0 / rad);
289 }
RichardsSumQuantity::zero
void zero()
sets _total = 0
Definition: RichardsSumQuantity.C:31
PeacemanBorehole::_zs
std::vector< Real > _zs
z points of borehole
Definition: PeacemanBorehole.h:91
PeacemanBorehole::parseNextLineReals
bool parseNextLineReals(std::ifstream &ifs, std::vector< Real > &myvec)
reads a space-separated line of floats from ifs and puts in myvec
Definition: PeacemanBorehole.C:158
PeacemanBorehole::_xs
std::vector< Real > _xs
x points of the borehole
Definition: PeacemanBorehole.h:85
PeacemanBorehole::_re_constant
const Real _re_constant
borehole constant
Definition: PeacemanBorehole.h:40
PeacemanBorehole::_rs
std::vector< Real > _rs
radii of the borehole
Definition: PeacemanBorehole.h:82
pow
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
Definition: ExpressionBuilder.h:673
PeacemanBorehole::wellConstant
Real wellConstant(const RealTensorValue &perm, const RealTensorValue &rot, const Real &half_len, const Elem *ele, const Real &rad)
Calculates Peaceman's form of the borehole well constant Z Chen, Y Zhang, Well flow models for variou...
Definition: PeacemanBorehole.C:194
validParams< PeacemanBorehole >
InputParameters validParams< PeacemanBorehole >()
Definition: PeacemanBorehole.C:17
PeacemanBorehole::_rot_matrix
std::vector< RealTensorValue > _rot_matrix
rotation matrix used in well_constant calculation
Definition: PeacemanBorehole.h:100
PeacemanBorehole::_bottom_point
Point _bottom_point
the bottom point of the borehole (where bottom_pressure is defined)
Definition: PeacemanBorehole.h:94
PeacemanBorehole.h
PeacemanBorehole::_half_seg_len
std::vector< Real > _half_seg_len
0.5*(length of polyline segments between points)
Definition: PeacemanBorehole.h:97
PeacemanBorehole::_ys
std::vector< Real > _ys
y points of the borehole
Definition: PeacemanBorehole.h:88
PeacemanBorehole::addPoints
virtual void addPoints()
Add Dirac Points to the borehole.
Definition: PeacemanBorehole.C:180
PeacemanBorehole::_well_constant
const Real _well_constant
well constant
Definition: PeacemanBorehole.h:43
PeacemanBorehole::_borehole_direction
const RealVectorValue _borehole_direction
borehole direction. Note this is only used if there is only one borehole point
Definition: PeacemanBorehole.h:49
PeacemanBorehole::_total_outflow_mass
RichardsSumQuantity & _total_outflow_mass
This is used to hold the total fluid flowing into the borehole Hence, it is positive for production w...
Definition: PeacemanBorehole.h:79
PeacemanBorehole::PeacemanBorehole
PeacemanBorehole(const InputParameters &parameters)
Creates a new PeacemanBorehole This reads the file containing the lines of the form radius x y z that...
Definition: PeacemanBorehole.C:80
RichardsSumQuantity
Sums into _total This is used, for instance, to record the total mass flowing into a borehole.
Definition: RichardsSumQuantity.h:26
PeacemanBorehole::_point_file
const std::string _point_file
File defining the geometry of the borehole.
Definition: PeacemanBorehole.h:56
PeacemanBorehole::_borehole_length
const Real _borehole_length
borehole length. Note this is only used if there is only one borehole point
Definition: PeacemanBorehole.h:46