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

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

#include <RichardsBorehole.h>

Inheritance diagram for RichardsBorehole:
[legend]

Public Member Functions

 RichardsBorehole (const InputParameters &parameters)
 Creates a new RichardsBorehole 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 Note: at March2014 this is never called since moose does not have this functionality. More...
 

Protected Member Functions

void prepareNodalValues ()
 calculates the nodal values of pressure, mobility, and derivatives thereof More...
 
Real jac (unsigned int wrt_num)
 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 bool _fully_upwind
 Whether to use full upwinding. More...
 
const RichardsVarNames_richards_name_UO
 Defines the richards variables in the simulation. More...
 
const unsigned int _num_p
 number of richards variables More...
 
const unsigned int _pvar
 The moose internal variable number of the richards variable of this Dirac Kernel. More...
 
const RichardsDensity_density_UO
 user object defining the density. Only used if _fully_upwind = true More...
 
const RichardsSeff_seff_UO
 user object defining the effective saturation. Only used if _fully_upwind = true More...
 
const RichardsRelPerm_relperm_UO
 user object defining the relative permeability. Only used if _fully_upwind = true More...
 
unsigned int _num_nodes
 number of nodes in this element. Only used if _fully_upwind = true More...
 
std::vector< Real > _mobility
 nodal values of mobility = density*relperm/viscosity These are used if _fully_upwind = true More...
 
std::vector< std::vector< Real > > _dmobility_dv
 d(_mobility)/d(variable_ph) (variable_ph is the variable for phase=ph) These are used in the jacobian calculations if _fully_upwind = true More...
 
const MaterialProperty< std::vector< Real > > & _pp
 fluid porepressure (or porepressures in case of multiphase) More...
 
const MaterialProperty< std::vector< std::vector< Real > > > & _dpp_dv
 d(porepressure_i)/d(variable_j) More...
 
const MaterialProperty< std::vector< Real > > & _viscosity
 fluid viscosity More...
 
const MaterialProperty< RealTensorValue > & _permeability
 material permeability More...
 
const MaterialProperty< std::vector< std::vector< Real > > > & _dseff_dv
 deriviatves of Seff wrt variables More...
 
const MaterialProperty< std::vector< Real > > & _rel_perm
 relative permeability More...
 
const MaterialProperty< std::vector< std::vector< Real > > > & _drel_perm_dv
 d(relperm_i)/d(variable_j) More...
 
const MaterialProperty< std::vector< Real > > & _density
 fluid density More...
 
const MaterialProperty< std::vector< std::vector< Real > > > & _ddensity_dv
 d(density_i)/d(variable_j) More...
 
std::vector< const VariableValue * > _ps_at_nodes
 Holds the values of pressures at all the nodes of the element Only used if _fully_upwind = true Eg: _ps_at_nodes[_pvar] is a pointer to this variable's nodal porepressure values So: (*_ps_at_nodes[_pvar])[i] = _var.nodalSln()[i] = porepressure of pressure-variable _pvar at node i. 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.

Definition at line 26 of file RichardsBorehole.h.

Constructor & Destructor Documentation

RichardsBorehole::RichardsBorehole ( const InputParameters &  parameters)

Creates a new RichardsBorehole 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 37 of file RichardsBorehole.C.

38  : PeacemanBorehole(parameters),
39  _fully_upwind(getParam<bool>("fully_upwind")),
40  _richards_name_UO(getUserObject<RichardsVarNames>("richardsVarNames_UO")),
43 
44  // in the following, getUserObjectByName returns a reference (an alias) to a RichardsBLAH user
45  // object, and the & turns it into a pointer
47  ? &getUserObjectByName<RichardsDensity>(
48  getParam<std::vector<UserObjectName>>("density_UO")[_pvar])
49  : NULL),
51  ? &getUserObjectByName<RichardsSeff>(
52  getParam<std::vector<UserObjectName>>("seff_UO")[_pvar])
53  : NULL),
55  ? &getUserObjectByName<RichardsRelPerm>(
56  getParam<std::vector<UserObjectName>>("relperm_UO")[_pvar])
57  : NULL),
58 
59  _num_nodes(0),
60  _mobility(0),
61  _dmobility_dv(0),
62  _pp(getMaterialProperty<std::vector<Real>>("porepressure")),
63  _dpp_dv(getMaterialProperty<std::vector<std::vector<Real>>>("dporepressure_dv")),
64  _viscosity(getMaterialProperty<std::vector<Real>>("viscosity")),
65  _permeability(getMaterialProperty<RealTensorValue>("permeability")),
66  _dseff_dv(getMaterialProperty<std::vector<std::vector<Real>>>("ds_eff_dv")),
67  _rel_perm(getMaterialProperty<std::vector<Real>>("rel_perm")),
68  _drel_perm_dv(getMaterialProperty<std::vector<std::vector<Real>>>("drel_perm_dv")),
69  _density(getMaterialProperty<std::vector<Real>>("density")),
70  _ddensity_dv(getMaterialProperty<std::vector<std::vector<Real>>>("ddensity_dv"))
71 {
72  _ps_at_nodes.resize(_num_p);
73  for (unsigned int pnum = 0; pnum < _num_p; ++pnum)
75 
76  // To correctly compute the Jacobian terms,
77  // tell MOOSE that this DiracKernel depends on all the Richards Vars
78  const std::vector<MooseVariable *> & coupled_vars = _richards_name_UO.getCoupledMooseVars();
79  for (unsigned int i = 0; i < coupled_vars.size(); i++)
80  addMooseVariableDependency(coupled_vars[i]);
81 }
const MaterialProperty< RealTensorValue > & _permeability
material permeability
const bool _fully_upwind
Whether to use full upwinding.
const RichardsDensity * _density_UO
user object defining the density. Only used if _fully_upwind = true
const MaterialProperty< std::vector< Real > > & _density
fluid density
PeacemanBorehole(const InputParameters &parameters)
Creates a new PeacemanBorehole This reads the file containing the lines of the form radius x y z that...
const MaterialProperty< std::vector< Real > > & _viscosity
fluid viscosity
const RichardsSeff * _seff_UO
user object defining the effective saturation. Only used if _fully_upwind = true
const MaterialProperty< std::vector< std::vector< Real > > > & _dseff_dv
deriviatves of Seff wrt variables
const VariableValue * nodal_var(unsigned int richards_var_num) const
The nodal variable values for the given richards_var_num To extract a the value of pressure variable ...
unsigned int _num_nodes
number of nodes in this element. Only used if _fully_upwind = true
std::vector< const VariableValue * > _ps_at_nodes
Holds the values of pressures at all the nodes of the element Only used if _fully_upwind = true Eg: _...
std::vector< Real > _mobility
nodal values of mobility = density*relperm/viscosity These are used if _fully_upwind = true ...
const MaterialProperty< std::vector< Real > > & _pp
fluid porepressure (or porepressures in case of multiphase)
const MaterialProperty< std::vector< Real > > & _rel_perm
relative permeability
const RichardsRelPerm * _relperm_UO
user object defining the relative permeability. Only used if _fully_upwind = true ...
unsigned int num_v() const
the number of porepressure variables
const MaterialProperty< std::vector< std::vector< Real > > > & _drel_perm_dv
d(relperm_i)/d(variable_j)
const unsigned int _pvar
The moose internal variable number of the richards variable of this Dirac Kernel. ...
std::vector< std::vector< Real > > _dmobility_dv
d(_mobility)/d(variable_ph) (variable_ph is the variable for phase=ph) These are used in the jacobian...
const MaterialProperty< std::vector< std::vector< Real > > > & _dpp_dv
d(porepressure_i)/d(variable_j)
const RichardsVarNames & _richards_name_UO
Defines the richards variables in the simulation.
const unsigned int _num_p
number of richards variables
const MaterialProperty< std::vector< std::vector< Real > > > & _ddensity_dv
d(density_i)/d(variable_j)
unsigned int richards_var_num(unsigned int moose_var_num) const
the richards variable number

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

Computes the Jacobian.

This just calls prepareNodalValues, if _fully_upwind then calls DiracKernel::computeJacobian

Definition at line 200 of file RichardsBorehole.C.

201 {
202  if (_fully_upwind)
204  DiracKernel::computeJacobian();
205 }
const bool _fully_upwind
Whether to use full upwinding.
void prepareNodalValues()
calculates the nodal values of pressure, mobility, and derivatives thereof
Real RichardsBorehole::computeQpJacobian ( )
virtual

Computes the diagonal part of the jacobian.

Definition at line 208 of file RichardsBorehole.C.

209 {
210  const Real character = _character.value(_t, _q_point[_qp]);
211  if (character == 0.0)
212  return 0.0;
213  return jac(_pvar);
214 }
Real jac(unsigned int wrt_num)
Calculates Jacobian.
const unsigned int _pvar
The moose internal variable number of the richards variable of this Dirac Kernel. ...
Function & _character
If positive then the borehole acts as a sink (producion well) for porepressure > borehole pressure...
Real RichardsBorehole::computeQpOffDiagJacobian ( unsigned int  jvar)
virtual

Computes the off-diagonal part of the jacobian Note: at March2014 this is never called since moose does not have this functionality.

Hence as of March2014 this has never been tested.

Definition at line 217 of file RichardsBorehole.C.

218 {
220  return 0.0;
221  const unsigned int dvar = _richards_name_UO.richards_var_num(jvar);
222  return jac(dvar);
223 }
Real jac(unsigned int wrt_num)
Calculates Jacobian.
bool not_richards_var(unsigned int moose_var_num) const
returns true if moose_var_num is not a richards var
const RichardsVarNames & _richards_name_UO
Defines the richards variables in the simulation.
unsigned int richards_var_num(unsigned int moose_var_num) const
the richards variable number
Real RichardsBorehole::computeQpResidual ( )
virtual

Computes the Qp residual.

Definition at line 132 of file RichardsBorehole.C.

133 {
134  const Real character = _character.value(_t, _q_point[_qp]);
135  if (character == 0.0)
136  return 0.0;
137 
138  const Real bh_pressure =
139  _p_bot +
140  _unit_weight *
141  (_q_point[_qp] -
142  _bottom_point); // really want to use _q_point instaed of _current_point, i think?!
143 
144  Real pp;
145  Real mob;
146  if (!_fully_upwind)
147  {
148  pp = _pp[_qp][_pvar];
149  mob = _rel_perm[_qp][_pvar] * _density[_qp][_pvar] / _viscosity[_qp][_pvar];
150  }
151  else
152  {
153  pp = (*_ps_at_nodes[_pvar])[_i];
154  mob = _mobility[_i];
155  }
156 
157  // Get the ID we initially assigned to this point
158  const unsigned current_dirac_ptid = currentPointCachedID();
159 
160  // If getting the ID failed, fall back to the old bodge!
161  // if (current_dirac_ptid == libMesh::invalid_uint)
162  // current_dirac_ptid = (_zs.size() > 2) ? 1 : 0;
163 
164  Real outflow(0.0); // this is the flow rate from porespace out of the system
165 
166  Real wc(0.0);
167  if (current_dirac_ptid > 0)
168  // contribution from half-segment "behind" this point (must have >1 point for
169  // current_dirac_ptid>0)
170  {
171  wc = wellConstant(_permeability[_qp],
172  _rot_matrix[current_dirac_ptid - 1],
173  _half_seg_len[current_dirac_ptid - 1],
174  _current_elem,
175  _rs[current_dirac_ptid]);
176  if ((character < 0.0 && pp < bh_pressure) || (character > 0.0 && pp > bh_pressure))
177  // injection, so outflow<0 || // production, so outflow>0
178  outflow += _test[_i][_qp] * std::abs(character) * wc * mob * (pp - bh_pressure);
179  }
180 
181  if (current_dirac_ptid + 1 < _zs.size() || _zs.size() == 1)
182  // contribution from half-segment "ahead of" this point, or we only have one point
183  {
184  wc = wellConstant(_permeability[_qp],
185  _rot_matrix[current_dirac_ptid],
186  _half_seg_len[current_dirac_ptid],
187  _current_elem,
188  _rs[current_dirac_ptid]);
189  if ((character < 0.0 && pp < bh_pressure) || (character > 0.0 && pp > bh_pressure))
190  // injection, so outflow<0 || // production, so outflow>0
191  outflow += _test[_i][_qp] * std::abs(character) * wc * mob * (pp - bh_pressure);
192  }
193 
195  outflow * _dt); // this is not thread safe, but DiracKernel's aren't currently threaded
196  return outflow;
197 }
const MaterialProperty< RealTensorValue > & _permeability
material permeability
const bool _fully_upwind
Whether to use full upwinding.
const MaterialProperty< std::vector< Real > > & _density
fluid density
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 MaterialProperty< std::vector< Real > > & _viscosity
fluid viscosity
std::vector< Real > _rs
radii of the borehole
std::vector< const VariableValue * > _ps_at_nodes
Holds the values of pressures at all the nodes of the element Only used if _fully_upwind = true Eg: _...
std::vector< Real > _mobility
nodal values of mobility = density*relperm/viscosity These are used if _fully_upwind = true ...
const MaterialProperty< std::vector< Real > > & _pp
fluid porepressure (or porepressures in case of multiphase)
const MaterialProperty< std::vector< Real > > & _rel_perm
relative permeability
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 unsigned int _pvar
The moose internal variable number of the richards variable of this Dirac Kernel. ...
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...
void RichardsBorehole::computeResidual ( )
virtual

Computes the residual.

This just calls prepareNodalValues, if _fully_upwind then calls DiracKernel::computeResidual

Definition at line 124 of file RichardsBorehole.C.

125 {
126  if (_fully_upwind)
128  DiracKernel::computeResidual();
129 }
const bool _fully_upwind
Whether to use full upwinding.
void prepareNodalValues()
calculates the nodal values of pressure, mobility, and derivatives thereof
Real RichardsBorehole::jac ( unsigned int  wrt_num)
protected

Calculates Jacobian.

Parameters
wrt_numdifferentiate the residual wrt this Richards variable

Definition at line 226 of file RichardsBorehole.C.

Referenced by computeQpJacobian(), and computeQpOffDiagJacobian().

227 {
228  const Real character = _character.value(_t, _q_point[_qp]);
229  if (character == 0.0)
230  return 0.0;
231 
232  const Real bh_pressure =
233  _p_bot +
234  _unit_weight *
235  (_q_point[_qp] -
236  _bottom_point); // really want to use _q_point instaed of _current_point, i think?!
237 
238  Real pp;
239  Real dpp_dv;
240  Real mob;
241  Real dmob_dv;
242  Real phi;
243  if (!_fully_upwind)
244  {
245  pp = _pp[_qp][_pvar];
246  dpp_dv = _dpp_dv[_qp][_pvar][wrt_num];
247  mob = _rel_perm[_qp][_pvar] * _density[_qp][_pvar] / _viscosity[_qp][_pvar];
248  dmob_dv = (_drel_perm_dv[_qp][_pvar][wrt_num] * _density[_qp][_pvar] +
249  _rel_perm[_qp][_pvar] * _ddensity_dv[_qp][_pvar][wrt_num]) /
250  _viscosity[_qp][_pvar];
251  phi = _phi[_j][_qp];
252  }
253  else
254  {
255  if (_i != _j)
256  return 0.0; // residual at node _i only depends on variables at that node
257  pp = (*_ps_at_nodes[_pvar])[_i];
258  dpp_dv =
259  (_pvar == wrt_num ? 1 : 0); // NOTE: i'm assuming that the variables are pressure variables
260  mob = _mobility[_i];
261  dmob_dv = _dmobility_dv[_i][wrt_num];
262  phi = 1;
263  }
264 
265  // Get the ID we initially assigned to this point
266  const unsigned current_dirac_ptid = currentPointCachedID();
267 
268  // If getting the ID failed, fall back to the old bodge!
269  // if (current_dirac_ptid == libMesh::invalid_uint)
270  // current_dirac_ptid = (_zs.size() > 2) ? 1 : 0;
271 
272  Real outflowp(0.0);
273 
274  Real wc(0.0);
275  if (current_dirac_ptid > 0)
276  // contribution from half-segment "behind" this point
277  {
278  wc = wellConstant(_permeability[_qp],
279  _rot_matrix[current_dirac_ptid - 1],
280  _half_seg_len[current_dirac_ptid - 1],
281  _current_elem,
282  _rs[current_dirac_ptid]);
283  if ((character < 0.0 && pp < bh_pressure) || (character > 0.0 && pp > bh_pressure))
284  // injection, so outflow<0 || // production, so outflow>0
285  outflowp += _test[_i][_qp] * std::abs(character) * wc *
286  (mob * phi * dpp_dv + dmob_dv * phi * (pp - bh_pressure));
287  }
288 
289  if (current_dirac_ptid < _zs.size() - 1 || _zs.size() == 1)
290  // contribution from half-segment "ahead of" this point
291  {
292  wc = wellConstant(_permeability[_qp],
293  _rot_matrix[current_dirac_ptid],
294  _half_seg_len[current_dirac_ptid],
295  _current_elem,
296  _rs[current_dirac_ptid]);
297  if ((character < 0.0 && pp < bh_pressure) || (character > 0.0 && pp > bh_pressure))
298  // injection, so outflow<0 || // production, so outflow>0
299  outflowp += _test[_i][_qp] * std::abs(character) * wc *
300  (mob * phi * dpp_dv + dmob_dv * phi * (pp - bh_pressure));
301  }
302 
303  return outflowp;
304 }
const MaterialProperty< RealTensorValue > & _permeability
material permeability
const bool _fully_upwind
Whether to use full upwinding.
const MaterialProperty< std::vector< Real > > & _density
fluid density
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 MaterialProperty< std::vector< Real > > & _viscosity
fluid viscosity
std::vector< Real > _rs
radii of the borehole
std::vector< const VariableValue * > _ps_at_nodes
Holds the values of pressures at all the nodes of the element Only used if _fully_upwind = true Eg: _...
std::vector< Real > _mobility
nodal values of mobility = density*relperm/viscosity These are used if _fully_upwind = true ...
const MaterialProperty< std::vector< Real > > & _pp
fluid porepressure (or porepressures in case of multiphase)
const MaterialProperty< std::vector< Real > > & _rel_perm
relative permeability
const RealVectorValue _unit_weight
unit weight of fluid in borehole (for calculating bottomhole pressure at each Dirac Point) ...
const MaterialProperty< std::vector< std::vector< Real > > > & _drel_perm_dv
d(relperm_i)/d(variable_j)
const Real _p_bot
bottomhole pressure of borehole
const unsigned int _pvar
The moose internal variable number of the richards variable of this Dirac Kernel. ...
std::vector< std::vector< Real > > _dmobility_dv
d(_mobility)/d(variable_ph) (variable_ph is the variable for phase=ph) These are used in the jacobian...
const MaterialProperty< std::vector< std::vector< Real > > > & _dpp_dv
d(porepressure_i)/d(variable_j)
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
const MaterialProperty< std::vector< std::vector< Real > > > & _ddensity_dv
d(density_i)/d(variable_j)
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 RichardsBorehole::prepareNodalValues ( )
protected

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

Definition at line 84 of file RichardsBorehole.C.

Referenced by computeJacobian(), and computeResidual().

85 {
86  // NOTE: i'm assuming that all the richards variables are pressure values
87 
88  _num_nodes = (*_ps_at_nodes[_pvar]).size();
89 
90  Real p;
91  Real density;
92  Real ddensity_dp;
93  Real seff;
94  std::vector<Real> dseff_dp;
95  Real relperm;
96  Real drelperm_ds;
97  _mobility.resize(_num_nodes);
98  _dmobility_dv.resize(_num_nodes);
99  dseff_dp.resize(_num_p);
100  for (unsigned int nodenum = 0; nodenum < _num_nodes; ++nodenum)
101  {
102  // retrieve and calculate basic things at the node
103  p = (*_ps_at_nodes[_pvar])[nodenum]; // pressure of fluid _pvar at node nodenum
104  density = _density_UO->density(p); // density of fluid _pvar at node nodenum
105  ddensity_dp = _density_UO->ddensity(p); // d(density)/dP
106  seff = _seff_UO->seff(_ps_at_nodes,
107  nodenum); // effective saturation of fluid _pvar at node nodenum
108  _seff_UO->dseff(
109  _ps_at_nodes, nodenum, dseff_dp); // d(seff)/d(P_ph), for ph = 0, ..., _num_p - 1
110  relperm = _relperm_UO->relperm(seff); // relative permeability of fluid _pvar at node nodenum
111  drelperm_ds = _relperm_UO->drelperm(seff); // d(relperm)/dseff
112 
113  // calculate the mobility and its derivatives wrt (variable_ph = porepressure_ph)
114  _mobility[nodenum] =
115  density * relperm / _viscosity[0][_pvar]; // assume viscosity is constant throughout element
116  _dmobility_dv[nodenum].resize(_num_p);
117  for (unsigned int ph = 0; ph < _num_p; ++ph)
118  _dmobility_dv[nodenum][ph] = density * drelperm_ds * dseff_dp[ph] / _viscosity[0][_pvar];
119  _dmobility_dv[nodenum][_pvar] += ddensity_dp * relperm / _viscosity[0][_pvar];
120  }
121 }
const RichardsDensity * _density_UO
user object defining the density. Only used if _fully_upwind = true
virtual Real drelperm(Real seff) const =0
derivative of relative permeability wrt effective saturation This must be over-ridden in your derived...
const MaterialProperty< std::vector< Real > > & _viscosity
fluid viscosity
const RichardsSeff * _seff_UO
user object defining the effective saturation. Only used if _fully_upwind = true
const std::string density
Definition: NS.h:15
unsigned int _num_nodes
number of nodes in this element. Only used if _fully_upwind = true
virtual void dseff(std::vector< const VariableValue * > p, unsigned int qp, std::vector< Real > &result) const =0
derivative(s) of effective saturation as a function of porepressure(s) at given quadpoint of the elem...
std::vector< const VariableValue * > _ps_at_nodes
Holds the values of pressures at all the nodes of the element Only used if _fully_upwind = true Eg: _...
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...
std::vector< Real > _mobility
nodal values of mobility = density*relperm/viscosity These are used if _fully_upwind = true ...
const RichardsRelPerm * _relperm_UO
user object defining the relative permeability. Only used if _fully_upwind = true ...
virtual Real ddensity(Real p) const =0
derivative of fluid density wrt porepressure This must be over-ridden in derived classes to provide a...
virtual Real relperm(Real seff) const =0
relative permeability as a function of effective saturation This must be over-ridden in your derived ...
const unsigned int _pvar
The moose internal variable number of the richards variable of this Dirac Kernel. ...
std::vector< std::vector< Real > > _dmobility_dv
d(_mobility)/d(variable_ph) (variable_ph is the variable for phase=ph) These are used in the jacobian...
const unsigned int _num_p
number of richards variables
virtual Real seff(std::vector< const VariableValue * > p, unsigned int qp) const =0
effective saturation as a function of porepressure(s) at given quadpoint of the element ...
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 Q2PBorehole::computeQpResidual(), computeQpResidual(), Q2PBorehole::jac(), and 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 Q2PBorehole::computeQpResidual(), computeQpResidual(), Q2PBorehole::jac(), 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 computeQpJacobian(), Q2PBorehole::computeQpResidual(), computeQpResidual(), Q2PBorehole::jac(), and jac().

const MaterialProperty<std::vector<std::vector<Real> > >& RichardsBorehole::_ddensity_dv
protected

d(density_i)/d(variable_j)

Definition at line 134 of file RichardsBorehole.h.

Referenced by jac().

const MaterialProperty<std::vector<Real> >& RichardsBorehole::_density
protected

fluid density

Definition at line 131 of file RichardsBorehole.h.

Referenced by computeQpResidual(), and jac().

const RichardsDensity* RichardsBorehole::_density_UO
protected

user object defining the density. Only used if _fully_upwind = true

Definition at line 86 of file RichardsBorehole.h.

Referenced by prepareNodalValues().

std::vector<std::vector<Real> > RichardsBorehole::_dmobility_dv
protected

d(_mobility)/d(variable_ph) (variable_ph is the variable for phase=ph) These are used in the jacobian calculations if _fully_upwind = true

Definition at line 107 of file RichardsBorehole.h.

Referenced by jac(), and prepareNodalValues().

const MaterialProperty<std::vector<std::vector<Real> > >& RichardsBorehole::_dpp_dv
protected

d(porepressure_i)/d(variable_j)

Definition at line 113 of file RichardsBorehole.h.

Referenced by jac().

const MaterialProperty<std::vector<std::vector<Real> > >& RichardsBorehole::_drel_perm_dv
protected

d(relperm_i)/d(variable_j)

Definition at line 128 of file RichardsBorehole.h.

Referenced by jac().

const MaterialProperty<std::vector<std::vector<Real> > >& RichardsBorehole::_dseff_dv
protected

deriviatves of Seff wrt variables

Definition at line 122 of file RichardsBorehole.h.

const bool RichardsBorehole::_fully_upwind
protected

Whether to use full upwinding.

Definition at line 74 of file RichardsBorehole.h.

Referenced by computeJacobian(), computeQpResidual(), computeResidual(), and jac().

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

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

nodal values of mobility = density*relperm/viscosity These are used if _fully_upwind = true

Definition at line 101 of file RichardsBorehole.h.

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

unsigned int RichardsBorehole::_num_nodes
protected

number of nodes in this element. Only used if _fully_upwind = true

Definition at line 95 of file RichardsBorehole.h.

Referenced by prepareNodalValues().

const unsigned int RichardsBorehole::_num_p
protected

number of richards variables

Definition at line 80 of file RichardsBorehole.h.

Referenced by prepareNodalValues(), and RichardsBorehole().

const Real PeacemanBorehole::_p_bot
protectedinherited

bottomhole pressure of borehole

Definition at line 68 of file PeacemanBorehole.h.

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

const MaterialProperty<RealTensorValue>& RichardsBorehole::_permeability
protected

material permeability

Definition at line 119 of file RichardsBorehole.h.

Referenced by computeQpResidual(), and jac().

const MaterialProperty<std::vector<Real> >& RichardsBorehole::_pp
protected

fluid porepressure (or porepressures in case of multiphase)

Definition at line 110 of file RichardsBorehole.h.

Referenced by computeQpResidual(), and jac().

std::vector<const VariableValue *> RichardsBorehole::_ps_at_nodes
protected

Holds the values of pressures at all the nodes of the element Only used if _fully_upwind = true Eg: _ps_at_nodes[_pvar] is a pointer to this variable's nodal porepressure values So: (*_ps_at_nodes[_pvar])[i] = _var.nodalSln()[i] = porepressure of pressure-variable _pvar at node i.

Definition at line 144 of file RichardsBorehole.h.

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

const unsigned int RichardsBorehole::_pvar
protected

The moose internal variable number of the richards variable of this Dirac Kernel.

Definition at line 83 of file RichardsBorehole.h.

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

const MaterialProperty<std::vector<Real> >& RichardsBorehole::_rel_perm
protected

relative permeability

Definition at line 125 of file RichardsBorehole.h.

Referenced by computeQpResidual(), and jac().

const RichardsRelPerm* RichardsBorehole::_relperm_UO
protected

user object defining the relative permeability. Only used if _fully_upwind = true

Definition at line 92 of file RichardsBorehole.h.

Referenced by prepareNodalValues().

const RichardsVarNames& RichardsBorehole::_richards_name_UO
protected

Defines the richards variables in the simulation.

Definition at line 77 of file RichardsBorehole.h.

Referenced by computeQpOffDiagJacobian(), and RichardsBorehole().

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

rotation matrix used in well_constant calculation

Definition at line 99 of file PeacemanBorehole.h.

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

std::vector<Real> PeacemanBorehole::_rs
protectedinherited
const RichardsSeff* RichardsBorehole::_seff_UO
protected

user object defining the effective saturation. Only used if _fully_upwind = true

Definition at line 89 of file RichardsBorehole.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(), Q2PBorehole::computeQpResidual(), 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 Q2PBorehole::computeQpResidual(), computeQpResidual(), Q2PBorehole::jac(), and jac().

const MaterialProperty<std::vector<Real> >& RichardsBorehole::_viscosity
protected

fluid viscosity

Definition at line 116 of file RichardsBorehole.h.

Referenced by computeQpResidual(), jac(), and 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: