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

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

#include <Q2PBorehole.h>

Inheritance diagram for Q2PBorehole:
[legend]

Public Member Functions

 Q2PBorehole (const InputParameters &parameters)
 Creates a new Q2PBorehole This sets all the variables, but also reads the file containing the lines of the form radius x y z that defines the borehole geometry. More...
 
virtual void computeResidual ()
 Computes the residual. More...
 
virtual Real computeQpResidual ()
 Computes the Qp residual. More...
 
virtual void computeJacobian ()
 Computes the Jacobian. More...
 
virtual Real computeQpJacobian ()
 Computes the diagonal part of the jacobian. More...
 
virtual Real computeQpOffDiagJacobian (unsigned int jvar)
 Computes the off-diagonal part of the jacobian. More...
 

Protected Member Functions

void prepareNodalValues ()
 calculates the nodal values of pressure, mobility, and derivatives thereof More...
 
Real jac (unsigned int jvar)
 Calculates Jacobian. More...
 
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

const RichardsDensity_density
 fluid density More...
 
const RichardsRelPerm_relperm
 fluid relative permeability More...
 
const VariableValue & _other_var_nodal
 the other variable in the 2-phase system (this is saturation if Variable=porepressure, and viceversa) More...
 
const unsigned int _other_var_num
 the variable number of the other variable More...
 
const bool _var_is_pp
 whether the Variable for this BC is porepressure or not More...
 
const Real _viscosity
 viscosity More...
 
const MaterialProperty< RealTensorValue > & _permeability
 permeability More...
 
unsigned int _num_nodes
 number of nodes in this element. More...
 
std::vector< Real > _pp
 nodal porepressure More...
 
std::vector< Real > _sat
 nodal saturation More...
 
std::vector< Real > _mobility
 nodal mobility More...
 
std::vector< Real > _dmobility_dp
 nodal d(mobility)/d(porepressure) More...
 
std::vector< Real > _dmobility_ds
 nodal d(mobility)/d(saturation) More...
 
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...
 

Detailed Description

Approximates a borehole by a sequence of Dirac Points.

This is for use by a Q2P model.

Definition at line 25 of file Q2PBorehole.h.

Constructor & Destructor Documentation

Q2PBorehole::Q2PBorehole ( const InputParameters &  parameters)

Creates a new Q2PBorehole This sets all the variables, but also 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 40 of file Q2PBorehole.C.

41  : PeacemanBorehole(parameters),
42  _density(getUserObject<RichardsDensity>("fluid_density")),
43  _relperm(getUserObject<RichardsRelPerm>("fluid_relperm")),
44  _other_var_nodal(coupledNodalValue("other_var")),
45  _other_var_num(coupled("other_var")),
46  _var_is_pp(getParam<bool>("var_is_porepressure")),
47  _viscosity(getParam<Real>("fluid_viscosity")),
48  _permeability(getMaterialProperty<RealTensorValue>("permeability")),
49  _num_nodes(0),
50  _pp(0),
51  _sat(0),
52  _mobility(0),
53  _dmobility_dp(0),
54  _dmobility_ds(0)
55 {
56 }
unsigned int _num_nodes
number of nodes in this element.
Definition: Q2PBorehole.h:91
const RichardsRelPerm & _relperm
fluid relative permeability
Definition: Q2PBorehole.h:73
PeacemanBorehole(const InputParameters &parameters)
Creates a new PeacemanBorehole This reads the file containing the lines of the form radius x y z that...
std::vector< Real > _dmobility_dp
nodal d(mobility)/d(porepressure)
Definition: Q2PBorehole.h:103
const unsigned int _other_var_num
the variable number of the other variable
Definition: Q2PBorehole.h:79
const VariableValue & _other_var_nodal
the other variable in the 2-phase system (this is saturation if Variable=porepressure, and viceversa)
Definition: Q2PBorehole.h:76
std::vector< Real > _mobility
nodal mobility
Definition: Q2PBorehole.h:100
const Real _viscosity
viscosity
Definition: Q2PBorehole.h:85
std::vector< Real > _sat
nodal saturation
Definition: Q2PBorehole.h:97
const RichardsDensity & _density
fluid density
Definition: Q2PBorehole.h:70
std::vector< Real > _dmobility_ds
nodal d(mobility)/d(saturation)
Definition: Q2PBorehole.h:106
const bool _var_is_pp
whether the Variable for this BC is porepressure or not
Definition: Q2PBorehole.h:82
const MaterialProperty< RealTensorValue > & _permeability
permeability
Definition: Q2PBorehole.h:88
std::vector< Real > _pp
nodal porepressure
Definition: Q2PBorehole.h:94

Member Function Documentation

void PeacemanBorehole::addPoints ( )
protectedvirtualinherited

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
void Q2PBorehole::computeJacobian ( )
virtual

Computes the Jacobian.

This just calls prepareNodalValues then calls DiracKernel::computeJacobian

Definition at line 167 of file Q2PBorehole.C.

168 {
170  DiracKernel::computeJacobian();
171 }
void prepareNodalValues()
calculates the nodal values of pressure, mobility, and derivatives thereof
Definition: Q2PBorehole.C:59
Real Q2PBorehole::computeQpJacobian ( )
virtual

Computes the diagonal part of the jacobian.

Definition at line 174 of file Q2PBorehole.C.

175 {
176  return jac(_var.number());
177 }
Real jac(unsigned int jvar)
Calculates Jacobian.
Definition: Q2PBorehole.C:188
Real Q2PBorehole::computeQpOffDiagJacobian ( unsigned int  jvar)
virtual

Computes the off-diagonal part of the jacobian.

Definition at line 180 of file Q2PBorehole.C.

181 {
182  if (jvar == _other_var_num || jvar == _var.number())
183  return jac(jvar);
184  return 0.0;
185 }
const unsigned int _other_var_num
the variable number of the other variable
Definition: Q2PBorehole.h:79
Real jac(unsigned int jvar)
Calculates Jacobian.
Definition: Q2PBorehole.C:188
Real Q2PBorehole::computeQpResidual ( )
virtual

Computes the Qp residual.

Definition at line 110 of file Q2PBorehole.C.

111 {
112  const Real character = _character.value(_t, _q_point[_qp]);
113  if (character == 0.0)
114  return 0.0;
115 
116  const Real bh_pressure =
117  _p_bot +
118  _unit_weight *
119  (_q_point[_qp] -
120  _bottom_point); // really want to use _q_point instaed of _current_point, i think?!
121 
122  // Get the ID we initially assigned to this point
123  const unsigned current_dirac_ptid = currentPointCachedID();
124 
125  // If getting the ID failed, fall back to the old bodge!
126  // if (current_dirac_ptid == libMesh::invalid_uint)
127  // current_dirac_ptid = (_zs.size() > 2) ? 1 : 0;
128 
129  Real outflow(0.0); // this is the flow rate from porespace out of the system
130 
131  Real wc(0.0);
132  if (current_dirac_ptid > 0)
133  // contribution from half-segment "behind" this point (must have >1 point for
134  // current_dirac_ptid>0)
135  {
136  wc = wellConstant(_permeability[0],
137  _rot_matrix[current_dirac_ptid - 1],
138  _half_seg_len[current_dirac_ptid - 1],
139  _current_elem,
140  _rs[current_dirac_ptid]);
141  if ((character < 0.0 && _pp[_i] < bh_pressure) || (character > 0.0 && _pp[_i] > bh_pressure))
142  // injection, so outflow<0 || // production, so outflow>0
143  outflow +=
144  _test[_i][_qp] * std::abs(character) * wc * _mobility[_i] * (_pp[_i] - bh_pressure);
145  }
146 
147  if (current_dirac_ptid + 1 < _zs.size() || _zs.size() == 1)
148  // contribution from half-segment "ahead of" this point, or we only have one point
149  {
150  wc = wellConstant(_permeability[0],
151  _rot_matrix[current_dirac_ptid],
152  _half_seg_len[current_dirac_ptid],
153  _current_elem,
154  _rs[current_dirac_ptid]);
155  if ((character < 0.0 && _pp[_i] < bh_pressure) || (character > 0.0 && _pp[_i] > bh_pressure))
156  // injection, so outflow<0 || // production, so outflow>0
157  outflow +=
158  _test[_i][_qp] * std::abs(character) * wc * _mobility[_i] * (_pp[_i] - bh_pressure);
159  }
160 
162  outflow * _dt); // this is not thread safe, but DiracKernel's aren't currently threaded
163  return outflow;
164 }
Real wellConstant(const RealTensorValue &perm, const RealTensorValue &rot, const Real &half_len, const Elem *ele, const Real &rad)
Calculates Peaceman&#39;s form of the borehole well constant Z Chen, Y Zhang, Well flow models for variou...
std::vector< Real > _rs
radii of the borehole
std::vector< Real > _mobility
nodal mobility
Definition: Q2PBorehole.h:100
const RealVectorValue _unit_weight
unit weight of fluid in borehole (for calculating bottomhole pressure at each Dirac Point) ...
void add(Real contrib)
adds contrib to _total
RichardsSumQuantity & _total_outflow_mass
This is used to hold the total fluid flowing into the borehole Hence, it is positive for production w...
const Real _p_bot
bottomhole pressure of borehole
const MaterialProperty< RealTensorValue > & _permeability
permeability
Definition: Q2PBorehole.h:88
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
std::vector< Real > _pp
nodal porepressure
Definition: Q2PBorehole.h:94
Function & _character
If positive then the borehole acts as a sink (producion well) for porepressure > borehole pressure...
void Q2PBorehole::computeResidual ( )
virtual

Computes the residual.

This just calls prepareNodalValues then calls DiracKernel::computeResidual

Definition at line 103 of file Q2PBorehole.C.

104 {
106  DiracKernel::computeResidual();
107 }
void prepareNodalValues()
calculates the nodal values of pressure, mobility, and derivatives thereof
Definition: Q2PBorehole.C:59
Real Q2PBorehole::jac ( unsigned int  jvar)
protected

Calculates Jacobian.

Parameters
jvardifferentiate the residual wrt this variable

Definition at line 188 of file Q2PBorehole.C.

Referenced by computeQpJacobian(), and computeQpOffDiagJacobian().

189 {
190  if (_i != _j)
191  return 0.0;
192 
193  const Real character = _character.value(_t, _q_point[_qp]);
194  if (character == 0.0)
195  return 0.0;
196 
197  const Real bh_pressure =
198  _p_bot +
199  _unit_weight *
200  (_q_point[_qp] -
201  _bottom_point); // really want to use _q_point instaed of _current_point, i think?!
202 
203  const Real phi = 1;
204 
205  // Get the ID we initially assigned to this point
206  const unsigned current_dirac_ptid = currentPointCachedID();
207 
208  // If getting the ID failed, fall back to the old bodge!
209  // if (current_dirac_ptid == libMesh::invalid_uint)
210  // current_dirac_ptid = (_zs.size() > 2) ? 1 : 0;
211 
212  Real outflowp(0.0);
213 
214  const bool deriv_wrt_pp =
215  (_var_is_pp && (jvar == _var.number())) || (!_var_is_pp && (jvar == _other_var_num));
216 
217  Real wc(0.0);
218  if (current_dirac_ptid > 0)
219  // contribution from half-segment "behind" this point
220  {
221  wc = wellConstant(_permeability[0],
222  _rot_matrix[current_dirac_ptid - 1],
223  _half_seg_len[current_dirac_ptid - 1],
224  _current_elem,
225  _rs[current_dirac_ptid]);
226  if ((character < 0.0 && _pp[_i] < bh_pressure) || (character > 0.0 && _pp[_i] > bh_pressure))
227  {
228  // injection, so outflow<0 || // production, so outflow>0
229  if (deriv_wrt_pp)
230  outflowp += std::abs(character) * wc *
231  (_mobility[_i] * phi + _dmobility_dp[_i] * phi * (_pp[_i] - bh_pressure));
232  else
233  outflowp += std::abs(character) * wc * _dmobility_ds[_i] * phi * (_pp[_i] - bh_pressure);
234  }
235  }
236 
237  if (current_dirac_ptid < _zs.size() - 1 || _zs.size() == 1)
238  // contribution from half-segment "ahead of" this point
239  {
240  wc = wellConstant(_permeability[0],
241  _rot_matrix[current_dirac_ptid],
242  _half_seg_len[current_dirac_ptid],
243  _current_elem,
244  _rs[current_dirac_ptid]);
245  if ((character < 0.0 && _pp[_i] < bh_pressure) || (character > 0.0 && _pp[_i] > bh_pressure))
246  {
247  // injection, so outflow<0 || // production, so outflow>0
248  if (deriv_wrt_pp)
249  outflowp += std::abs(character) * wc *
250  (_mobility[_i] * phi + _dmobility_dp[_i] * phi * (_pp[_i] - bh_pressure));
251  else
252  outflowp += std::abs(character) * wc * _dmobility_ds[_i] * phi * (_pp[_i] - bh_pressure);
253  }
254  }
255 
256  return _test[_i][_qp] * outflowp;
257 }
std::vector< Real > _dmobility_dp
nodal d(mobility)/d(porepressure)
Definition: Q2PBorehole.h:103
Real wellConstant(const RealTensorValue &perm, const RealTensorValue &rot, const Real &half_len, const Elem *ele, const Real &rad)
Calculates Peaceman&#39;s form of the borehole well constant Z Chen, Y Zhang, Well flow models for variou...
const unsigned int _other_var_num
the variable number of the other variable
Definition: Q2PBorehole.h:79
std::vector< Real > _rs
radii of the borehole
std::vector< Real > _mobility
nodal mobility
Definition: Q2PBorehole.h:100
const RealVectorValue _unit_weight
unit weight of fluid in borehole (for calculating bottomhole pressure at each Dirac Point) ...
const Real _p_bot
bottomhole pressure of borehole
std::vector< Real > _dmobility_ds
nodal d(mobility)/d(saturation)
Definition: Q2PBorehole.h:106
const bool _var_is_pp
whether the Variable for this BC is porepressure or not
Definition: Q2PBorehole.h:82
const MaterialProperty< RealTensorValue > & _permeability
permeability
Definition: Q2PBorehole.h:88
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
std::vector< Real > _pp
nodal porepressure
Definition: Q2PBorehole.h:94
Function & _character
If positive then the borehole acts as a sink (producion well) for porepressure > borehole pressure...
bool PeacemanBorehole::parseNextLineReals ( std::ifstream &  ifs,
std::vector< Real > &  myvec 
)
protectedinherited

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

Definition at line 156 of file PeacemanBorehole.C.

Referenced by PeacemanBorehole::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 }
void Q2PBorehole::prepareNodalValues ( )
protected

calculates the nodal values of pressure, mobility, and derivatives thereof

Definition at line 59 of file Q2PBorehole.C.

Referenced by computeJacobian(), and computeResidual().

60 {
62 
63  // set _pp and _sat variables
64  _pp.resize(_num_nodes);
65  _sat.resize(_num_nodes);
66  if (_var_is_pp)
67  {
68  for (unsigned int nodenum = 0; nodenum < _num_nodes; ++nodenum)
69  {
70  _pp[nodenum] = _var.nodalSln()[nodenum];
71  _sat[nodenum] = _other_var_nodal[nodenum];
72  }
73  }
74  else
75  {
76  for (unsigned int nodenum = 0; nodenum < _num_nodes; ++nodenum)
77  {
78  _pp[nodenum] = _other_var_nodal[nodenum];
79  _sat[nodenum] = _var.nodalSln()[nodenum];
80  }
81  }
82 
83  Real density;
84  Real ddensity_dp;
85  Real relperm;
86  Real drelperm_ds;
87  _mobility.resize(_num_nodes);
88  _dmobility_dp.resize(_num_nodes);
89  _dmobility_ds.resize(_num_nodes);
90  for (unsigned int nodenum = 0; nodenum < _num_nodes; ++nodenum)
91  {
92  density = _density.density(_pp[nodenum]);
93  ddensity_dp = _density.ddensity(_pp[nodenum]);
94  relperm = _relperm.relperm(_sat[nodenum]);
95  drelperm_ds = _relperm.drelperm(_sat[nodenum]);
96  _mobility[nodenum] = density * relperm / _viscosity;
97  _dmobility_dp[nodenum] = ddensity_dp * relperm / _viscosity;
98  _dmobility_ds[nodenum] = density * drelperm_ds / _viscosity;
99  }
100 }
unsigned int _num_nodes
number of nodes in this element.
Definition: Q2PBorehole.h:91
const RichardsRelPerm & _relperm
fluid relative permeability
Definition: Q2PBorehole.h:73
virtual Real drelperm(Real seff) const =0
derivative of relative permeability wrt effective saturation This must be over-ridden in your derived...
std::vector< Real > _dmobility_dp
nodal d(mobility)/d(porepressure)
Definition: Q2PBorehole.h:103
const std::string density
Definition: NS.h:15
const VariableValue & _other_var_nodal
the other variable in the 2-phase system (this is saturation if Variable=porepressure, and viceversa)
Definition: Q2PBorehole.h:76
std::vector< Real > _mobility
nodal mobility
Definition: Q2PBorehole.h:100
virtual Real density(Real p) const =0
fluid density as a function of porepressure This must be over-ridden in derived classes to provide an...
const Real _viscosity
viscosity
Definition: Q2PBorehole.h:85
std::vector< Real > _sat
nodal saturation
Definition: Q2PBorehole.h:97
const RichardsDensity & _density
fluid density
Definition: Q2PBorehole.h:70
virtual Real ddensity(Real p) const =0
derivative of fluid density wrt porepressure This must be over-ridden in derived classes to provide a...
std::vector< Real > _dmobility_ds
nodal d(mobility)/d(saturation)
Definition: Q2PBorehole.h:106
const bool _var_is_pp
whether the Variable for this BC is porepressure or not
Definition: Q2PBorehole.h:82
virtual Real relperm(Real seff) const =0
relative permeability as a function of effective saturation This must be over-ridden in your derived ...
std::vector< Real > _pp
nodal porepressure
Definition: Q2PBorehole.h:94
Real PeacemanBorehole::wellConstant ( const RealTensorValue &  perm,
const RealTensorValue &  rot,
const Real &  half_len,
const Elem *  ele,
const Real &  rad 
)
protectedinherited

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 computeQpResidual(), RichardsBorehole::computeQpResidual(), 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

Point PeacemanBorehole::_bottom_point
protectedinherited

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

Definition at line 93 of file PeacemanBorehole.h.

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

Function& PeacemanBorehole::_character
protectedinherited

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(), computeQpResidual(), RichardsBorehole::computeQpResidual(), jac(), and RichardsBorehole::jac().

const RichardsDensity& Q2PBorehole::_density
protected

fluid density

Definition at line 70 of file Q2PBorehole.h.

Referenced by prepareNodalValues().

std::vector<Real> Q2PBorehole::_dmobility_dp
protected

nodal d(mobility)/d(porepressure)

Definition at line 103 of file Q2PBorehole.h.

Referenced by jac(), and prepareNodalValues().

std::vector<Real> Q2PBorehole::_dmobility_ds
protected

nodal d(mobility)/d(saturation)

Definition at line 106 of file Q2PBorehole.h.

Referenced by jac(), and prepareNodalValues().

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

0.5*(length of polyline segments between points)

Definition at line 96 of file PeacemanBorehole.h.

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

std::vector<Real> Q2PBorehole::_mobility
protected

nodal mobility

Definition at line 100 of file Q2PBorehole.h.

Referenced by computeQpResidual(), jac(), and prepareNodalValues().

unsigned int Q2PBorehole::_num_nodes
protected

number of nodes in this element.

Definition at line 91 of file Q2PBorehole.h.

Referenced by prepareNodalValues().

const VariableValue& Q2PBorehole::_other_var_nodal
protected

the other variable in the 2-phase system (this is saturation if Variable=porepressure, and viceversa)

Definition at line 76 of file Q2PBorehole.h.

Referenced by prepareNodalValues().

const unsigned int Q2PBorehole::_other_var_num
protected

the variable number of the other variable

Definition at line 79 of file Q2PBorehole.h.

Referenced by computeQpOffDiagJacobian(), and jac().

const Real PeacemanBorehole::_p_bot
protectedinherited

bottomhole pressure of borehole

Definition at line 68 of file PeacemanBorehole.h.

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

const MaterialProperty<RealTensorValue>& Q2PBorehole::_permeability
protected

permeability

Definition at line 88 of file Q2PBorehole.h.

Referenced by computeQpResidual(), and jac().

std::vector<Real> Q2PBorehole::_pp
protected

nodal porepressure

Definition at line 94 of file Q2PBorehole.h.

Referenced by computeQpResidual(), jac(), and prepareNodalValues().

const RichardsRelPerm& Q2PBorehole::_relperm
protected

fluid relative permeability

Definition at line 73 of file Q2PBorehole.h.

Referenced by prepareNodalValues().

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

rotation matrix used in well_constant calculation

Definition at line 99 of file PeacemanBorehole.h.

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

std::vector<Real> PeacemanBorehole::_rs
protectedinherited
std::vector<Real> Q2PBorehole::_sat
protected

nodal saturation

Definition at line 97 of file Q2PBorehole.h.

Referenced by prepareNodalValues().

RichardsSumQuantity& PeacemanBorehole::_total_outflow_mass
protectedinherited

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 PeacemanBorehole::addPoints(), computeQpResidual(), RichardsBorehole::computeQpResidual(), and PeacemanBorehole::PeacemanBorehole().

const RealVectorValue PeacemanBorehole::_unit_weight
protectedinherited

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

Definition at line 71 of file PeacemanBorehole.h.

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

const bool Q2PBorehole::_var_is_pp
protected

whether the Variable for this BC is porepressure or not

Definition at line 82 of file Q2PBorehole.h.

Referenced by jac(), and prepareNodalValues().

const Real Q2PBorehole::_viscosity
protected

viscosity

Definition at line 85 of file Q2PBorehole.h.

Referenced by prepareNodalValues().

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

x points of the borehole

Definition at line 84 of file PeacemanBorehole.h.

Referenced by PeacemanBorehole::addPoints(), and PeacemanBorehole::PeacemanBorehole().

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

y points of the borehole

Definition at line 87 of file PeacemanBorehole.h.

Referenced by PeacemanBorehole::addPoints(), and PeacemanBorehole::PeacemanBorehole().

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

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