www.mooseframework.org
Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
PeacemanBorehole Class Reference

Approximates a borehole by a sequence of Dirac Points. More...

#include <PeacemanBorehole.h>

Inheritance diagram for PeacemanBorehole:
[legend]

Public Member Functions

 PeacemanBorehole (const InputParameters &parameters)
 Creates a new PeacemanBorehole This reads the file containing the lines of the form radius x y z that defines the borehole geometry. More...
 

Protected Member Functions

virtual void addPoints ()
 Add Dirac Points to the borehole. More...
 
bool parseNextLineReals (std::ifstream &ifs, std::vector< Real > &myvec)
 reads a space-separated line of floats from ifs and puts in myvec More...
 
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 various numerical methods, Int J Num Analysis and Modeling, 3 (2008) 375-388. More...
 

Protected Attributes

Function & _character
 If positive then the borehole acts as a sink (producion well) for porepressure > borehole pressure, and does nothing otherwise If negative then the borehole acts as a source (injection well) for porepressure < borehole pressure, and does nothing otherwise The flow rate to/from the borehole is multiplied by |character|, so usually character = +/- 1. More...
 
const Real _p_bot
 bottomhole pressure of borehole More...
 
const RealVectorValue _unit_weight
 unit weight of fluid in borehole (for calculating bottomhole pressure at each Dirac Point) More...
 
RichardsSumQuantity_total_outflow_mass
 This is used to hold the total fluid flowing into the borehole Hence, it is positive for production wells where fluid is flowing from porespace into the borehole and removed from the model. More...
 
std::vector< Real > _rs
 radii of the borehole More...
 
std::vector< Real > _xs
 x points of the borehole More...
 
std::vector< Real > _ys
 y points of the borehole More...
 
std::vector< Real > _zs
 z points of borehole More...
 
Point _bottom_point
 the bottom point of the borehole (where bottom_pressure is defined) More...
 
std::vector< Real > _half_seg_len
 0.5*(length of polyline segments between points) More...
 
std::vector< RealTensorValue > _rot_matrix
 rotation matrix used in well_constant calculation More...
 

Private Attributes

const Real _re_constant
 borehole constant More...
 
const Real _well_constant
 well constant More...
 
const Real _borehole_length
 borehole length. Note this is only used if there is only one borehole point More...
 
const RealVectorValue _borehole_direction
 borehole direction. Note this is only used if there is only one borehole point More...
 
const std::string _point_file
 File defining the geometry of the borehole. More...
 

Detailed Description

Approximates a borehole by a sequence of Dirac Points.

Definition at line 24 of file PeacemanBorehole.h.

Constructor & Destructor Documentation

PeacemanBorehole::PeacemanBorehole ( const InputParameters &  parameters)

Creates a new PeacemanBorehole This reads the file containing the lines of the form radius x y z that defines the borehole geometry.

It also calculates segment-lengths and rotation matrices needed for computing the borehole well constant

Definition at line 78 of file PeacemanBorehole.C.

79  : DiracKernel(parameters),
80  _re_constant(getParam<Real>("re_constant")),
81  _well_constant(getParam<Real>("well_constant")),
82  _borehole_length(getParam<Real>("borehole_length")),
83  _borehole_direction(getParam<RealVectorValue>("borehole_direction")),
84  _point_file(getParam<std::string>("point_file")),
85  _character(getFunction("character")),
86  _p_bot(getParam<Real>("bottom_pressure")),
87  _unit_weight(getParam<RealVectorValue>("unit_weight")),
89  const_cast<RichardsSumQuantity &>(getUserObject<RichardsSumQuantity>("SumQuantityUO")))
90 {
91  // zero the outflow mass
93 
94  // open file
95  std::ifstream file(_point_file.c_str());
96  if (!file.good())
97  mooseError("Error opening file '" + _point_file + "' from a Peaceman-type Borehole.");
98 
99  // construct the arrays of radius, x, y and z
100  std::vector<Real> scratch;
101  while (parseNextLineReals(file, scratch))
102  {
103  if (scratch.size() >= 2)
104  {
105  _rs.push_back(scratch[0]);
106  _xs.push_back(scratch[1]);
107  if (scratch.size() >= 3)
108  _ys.push_back(scratch[2]);
109  else
110  _ys.push_back(0.0);
111  if (scratch.size() >= 4)
112  _zs.push_back(scratch[3]);
113  else
114  _zs.push_back(0.0);
115  }
116  }
117 
118  file.close();
119 
120  const int num_pts = _zs.size();
121  _bottom_point(0) = _xs[num_pts - 1];
122  _bottom_point(1) = _ys[num_pts - 1];
123  _bottom_point(2) = _zs[num_pts - 1];
124 
125  // construct the line-segment lengths between each point
126  _half_seg_len.resize(std::max(num_pts - 1, 1));
127  for (unsigned int i = 0; i + 1 < _xs.size(); ++i)
128  {
129  _half_seg_len[i] =
130  0.5 * std::sqrt(std::pow(_xs[i + 1] - _xs[i], 2) + std::pow(_ys[i + 1] - _ys[i], 2) +
131  std::pow(_zs[i + 1] - _zs[i], 2));
132  if (_half_seg_len[i] == 0)
133  mooseError("Peaceman-type borehole has a zero-segment length at (x,y,z) = ",
134  _xs[i],
135  " ",
136  _ys[i],
137  " ",
138  _zs[i],
139  "\n");
140  }
141  if (num_pts == 1)
143 
144  // construct the rotation matrix needed to rotate the permeability
145  _rot_matrix.resize(std::max(num_pts - 1, 1));
146  for (unsigned int i = 0; i + 1 < _xs.size(); ++i)
147  {
148  const RealVectorValue v2(_xs[i + 1] - _xs[i], _ys[i + 1] - _ys[i], _zs[i + 1] - _zs[i]);
149  _rot_matrix[i] = RotationMatrix::rotVecToZ(v2);
150  }
151  if (num_pts == 1)
152  _rot_matrix[0] = RotationMatrix::rotVecToZ(_borehole_direction);
153 }
const Real _well_constant
well constant
const RealVectorValue _borehole_direction
borehole direction. Note this is only used if there is only one borehole point
const Real _re_constant
borehole constant
std::vector< Real > _rs
radii of the borehole
bool parseNextLineReals(std::ifstream &ifs, std::vector< Real > &myvec)
reads a space-separated line of floats from ifs and puts in myvec
std::vector< Real > _xs
x points of the borehole
std::vector< Real > _ys
y points of the borehole
void zero()
sets _total = 0
const RealVectorValue _unit_weight
unit weight of fluid in borehole (for calculating bottomhole pressure at each Dirac Point) ...
RichardsSumQuantity & _total_outflow_mass
This is used to hold the total fluid flowing into the borehole Hence, it is positive for production w...
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
const Real _p_bot
bottomhole pressure of borehole
const Real _borehole_length
borehole length. Note this is only used if there is only one borehole point
const std::string _point_file
File defining the geometry of the borehole.
Point _bottom_point
the bottom point of the borehole (where bottom_pressure is defined)
std::vector< RealTensorValue > _rot_matrix
rotation matrix used in well_constant calculation
std::vector< Real > _half_seg_len
0.5*(length of polyline segments between points)
std::vector< Real > _zs
z points of borehole
Function & _character
If positive then the borehole acts as a sink (producion well) for porepressure > borehole pressure...

Member Function Documentation

void PeacemanBorehole::addPoints ( )
protectedvirtual

Add Dirac Points to the borehole.

Definition at line 178 of file PeacemanBorehole.C.

179 {
180  // This function gets called just before the DiracKernel is evaluated
181  // so this is a handy place to zero this out.
183 
184  // Add point using the unique ID "i", let the DiracKernel take
185  // care of the caching. This should be fast after the first call,
186  // as long as the points don't move around.
187  for (unsigned int i = 0; i < _zs.size(); i++)
188  addPoint(Point(_xs[i], _ys[i], _zs[i]), i);
189 }
std::vector< Real > _xs
x points of the borehole
std::vector< Real > _ys
y points of the borehole
void zero()
sets _total = 0
RichardsSumQuantity & _total_outflow_mass
This is used to hold the total fluid flowing into the borehole Hence, it is positive for production w...
std::vector< Real > _zs
z points of borehole
bool PeacemanBorehole::parseNextLineReals ( std::ifstream &  ifs,
std::vector< Real > &  myvec 
)
protected

reads a space-separated line of floats from ifs and puts in myvec

Definition at line 156 of file PeacemanBorehole.C.

Referenced by PeacemanBorehole().

158 {
159  std::string line;
160  myvec.clear();
161  bool gotline(false);
162  if (getline(ifs, line))
163  {
164  gotline = true;
165 
166  // Harvest floats separated by whitespace
167  std::istringstream iss(line);
168  Real f;
169  while (iss >> f)
170  {
171  myvec.push_back(f);
172  }
173  }
174  return gotline;
175 }
Real PeacemanBorehole::wellConstant ( const RealTensorValue &  perm,
const RealTensorValue &  rot,
const Real &  half_len,
const Elem *  ele,
const Real &  rad 
)
protected

Calculates Peaceman's form of the borehole well constant Z Chen, Y Zhang, Well flow models for various numerical methods, Int J Num Analysis and Modeling, 3 (2008) 375-388.

Definition at line 192 of file PeacemanBorehole.C.

Referenced by Q2PBorehole::computeQpResidual(), RichardsBorehole::computeQpResidual(), Q2PBorehole::jac(), and RichardsBorehole::jac().

198 {
199  if (_well_constant > 0)
200  return _well_constant;
201 
202  // rot_perm has its "2" component lying along the half segment
203  // we want to determine the eigenvectors of rot(0:1, 0:1), since, when
204  // rotated back to the original frame we will determine the element
205  // lengths along these directions
206  const RealTensorValue rot_perm = (rot * perm) * rot.transpose();
207  const Real trace2D = rot_perm(0, 0) + rot_perm(1, 1);
208  const Real det2D = rot_perm(0, 0) * rot_perm(1, 1) - rot_perm(0, 1) * rot_perm(1, 0);
209  const Real sq = std::sqrt(std::max(0.25 * trace2D * trace2D - det2D,
210  0.0)); // the std::max accounts for wierdo precision loss
211  const Real eig_val1 = 0.5 * trace2D + sq;
212  const Real eig_val2 = 0.5 * trace2D - sq;
213  RealVectorValue eig_vec1, eig_vec2;
214  if (sq > std::abs(trace2D) * 1E-7) // matrix is not a multiple of the identity (1E-7 accounts for
215  // precision in a crude way)
216  {
217  if (rot_perm(1, 0) != 0)
218  {
219  eig_vec1(0) = eig_val1 - rot_perm(1, 1);
220  eig_vec1(1) = rot_perm(1, 0);
221  eig_vec2(0) = eig_val2 - rot_perm(1, 1);
222  eig_vec2(1) = rot_perm(1, 0);
223  }
224  else if (rot_perm(0, 1) != 0)
225  {
226  eig_vec1(0) = rot_perm(0, 1);
227  eig_vec1(1) = eig_val1 - rot_perm(0, 0);
228  eig_vec2(0) = rot_perm(0, 1);
229  eig_vec2(1) = eig_val2 - rot_perm(0, 0);
230  }
231  else // off diagonal terms are both zero
232  {
233  eig_vec1(0) = 1;
234  eig_vec2(1) = 1;
235  }
236  }
237  else // matrix is basically a multiple of the identity
238  {
239  eig_vec1(0) = 1;
240  eig_vec2(1) = 1;
241  }
242 
243  // finally, rotate these to original frame and normalise
244  eig_vec1 = rot.transpose() * eig_vec1;
245  eig_vec1 /= std::sqrt(eig_vec1 * eig_vec1);
246  eig_vec2 = rot.transpose() * eig_vec2;
247  eig_vec2 /= std::sqrt(eig_vec2 * eig_vec2);
248 
249  // find the "length" of the element in these directions
250  // TODO - probably better to use variance than max&min
251  Real max1 = eig_vec1 * ele->point(0);
252  Real max2 = eig_vec2 * ele->point(0);
253  Real min1 = max1;
254  Real min2 = max2;
255  Real proj;
256  for (unsigned int i = 1; i < ele->n_nodes(); i++)
257  {
258  proj = eig_vec1 * ele->point(i);
259  max1 = (max1 < proj) ? proj : max1;
260  min1 = (min1 < proj) ? min1 : proj;
261 
262  proj = eig_vec2 * ele->point(i);
263  max2 = (max2 < proj) ? proj : max2;
264  min2 = (min2 < proj) ? min2 : proj;
265  }
266  const Real ll1 = max1 - min1;
267  const Real ll2 = max2 - min2;
268 
269  const Real r0 = _re_constant * std::sqrt(std::sqrt(eig_val1 / eig_val2) * std::pow(ll2, 2) +
270  std::sqrt(eig_val2 / eig_val1) * std::pow(ll1, 2)) /
271  (std::pow(eig_val1 / eig_val2, 0.25) + std::pow(eig_val2 / eig_val1, 0.25));
272 
273  const Real effective_perm = std::sqrt(det2D);
274 
275  const Real halfPi = acos(0.0);
276 
277  if (r0 <= rad)
278  mooseError("The effective element size (about 0.2-times-true-ele-size) for an element "
279  "containing a Peaceman-type borehole must be (much) larger than the borehole radius "
280  "for the Peaceman formulation to be correct. Your element has effective size ",
281  r0,
282  " and the borehole radius is ",
283  rad,
284  "\n");
285 
286  return 4 * halfPi * effective_perm * half_len / std::log(r0 / rad);
287 }
const Real _well_constant
well constant
const Real _re_constant
borehole constant
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)

Member Data Documentation

const RealVectorValue PeacemanBorehole::_borehole_direction
private

borehole direction. Note this is only used if there is only one borehole point

Definition at line 48 of file PeacemanBorehole.h.

Referenced by PeacemanBorehole().

const Real PeacemanBorehole::_borehole_length
private

borehole length. Note this is only used if there is only one borehole point

Definition at line 45 of file PeacemanBorehole.h.

Referenced by PeacemanBorehole().

Point PeacemanBorehole::_bottom_point
protected

the bottom point of the borehole (where bottom_pressure is defined)

Definition at line 93 of file PeacemanBorehole.h.

Referenced by Q2PBorehole::computeQpResidual(), RichardsBorehole::computeQpResidual(), Q2PBorehole::jac(), RichardsBorehole::jac(), and PeacemanBorehole().

Function& PeacemanBorehole::_character
protected

If positive then the borehole acts as a sink (producion well) for porepressure > borehole pressure, and does nothing otherwise If negative then the borehole acts as a source (injection well) for porepressure < borehole pressure, and does nothing otherwise The flow rate to/from the borehole is multiplied by |character|, so usually character = +/- 1.

Definition at line 65 of file PeacemanBorehole.h.

Referenced by RichardsBorehole::computeQpJacobian(), Q2PBorehole::computeQpResidual(), RichardsBorehole::computeQpResidual(), Q2PBorehole::jac(), and RichardsBorehole::jac().

std::vector<Real> PeacemanBorehole::_half_seg_len
protected

0.5*(length of polyline segments between points)

Definition at line 96 of file PeacemanBorehole.h.

Referenced by Q2PBorehole::computeQpResidual(), RichardsBorehole::computeQpResidual(), Q2PBorehole::jac(), RichardsBorehole::jac(), and PeacemanBorehole().

const Real PeacemanBorehole::_p_bot
protected

bottomhole pressure of borehole

Definition at line 68 of file PeacemanBorehole.h.

Referenced by Q2PBorehole::computeQpResidual(), RichardsBorehole::computeQpResidual(), Q2PBorehole::jac(), and RichardsBorehole::jac().

const std::string PeacemanBorehole::_point_file
private

File defining the geometry of the borehole.

Each row has format radius x y z and the list of such points defines a polyline that is the borehole

Definition at line 55 of file PeacemanBorehole.h.

Referenced by PeacemanBorehole().

const Real PeacemanBorehole::_re_constant
private

borehole constant

Definition at line 39 of file PeacemanBorehole.h.

Referenced by wellConstant().

std::vector<RealTensorValue> PeacemanBorehole::_rot_matrix
protected

rotation matrix used in well_constant calculation

Definition at line 99 of file PeacemanBorehole.h.

Referenced by Q2PBorehole::computeQpResidual(), RichardsBorehole::computeQpResidual(), Q2PBorehole::jac(), RichardsBorehole::jac(), and PeacemanBorehole().

std::vector<Real> PeacemanBorehole::_rs
protected
RichardsSumQuantity& PeacemanBorehole::_total_outflow_mass
protected

This is used to hold the total fluid flowing into the borehole Hence, it is positive for production wells where fluid is flowing from porespace into the borehole and removed from the model.

Definition at line 78 of file PeacemanBorehole.h.

Referenced by addPoints(), Q2PBorehole::computeQpResidual(), RichardsBorehole::computeQpResidual(), and PeacemanBorehole().

const RealVectorValue PeacemanBorehole::_unit_weight
protected

unit weight of fluid in borehole (for calculating bottomhole pressure at each Dirac Point)

Definition at line 71 of file PeacemanBorehole.h.

Referenced by Q2PBorehole::computeQpResidual(), RichardsBorehole::computeQpResidual(), Q2PBorehole::jac(), and RichardsBorehole::jac().

const Real PeacemanBorehole::_well_constant
private

well constant

Definition at line 42 of file PeacemanBorehole.h.

Referenced by wellConstant().

std::vector<Real> PeacemanBorehole::_xs
protected

x points of the borehole

Definition at line 84 of file PeacemanBorehole.h.

Referenced by addPoints(), and PeacemanBorehole().

std::vector<Real> PeacemanBorehole::_ys
protected

y points of the borehole

Definition at line 87 of file PeacemanBorehole.h.

Referenced by addPoints(), and PeacemanBorehole().

std::vector<Real> PeacemanBorehole::_zs
protected

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