Go to the documentation of this file.
11 #include "RotationMatrix.h"
20 params.addRequiredParam<UserObjectName>(
22 "A RichardsDensity UserObject that defines the fluid density as a function of pressure.");
23 params.addRequiredParam<UserObjectName>(
25 "A RichardsRelPerm UserObject (eg RichardsRelPermPower) that defines the "
26 "fluid relative permeability as a function of the saturation Variable.");
27 params.addRequiredCoupledVar(
"other_var",
28 "The other variable in the 2-phase system. If "
29 "Variable=porepressure, the other_var=saturation, and "
31 params.addRequiredParam<
bool>(
"var_is_porepressure",
32 "This flag is needed to correctly calculate the Jacobian entries. "
33 "If set to true, this Sink will extract fluid from the phase with "
34 "porepressure as its Variable (usually the liquid phase). If set "
35 "to false, this Sink will extract fluid from the phase with "
36 "saturation as its variable (usually the gas phase)");
37 params.addRequiredParam<Real>(
"fluid_viscosity",
"The fluid dynamic viscosity");
38 params.addClassDescription(
"Approximates a borehole in the mesh with given bottomhole pressure, "
39 "and radii using a number of point sinks whose positions are read "
40 "from a file. This DiracKernel is for use by Q2P models");
48 _other_var_nodal(coupledDofValues(
"other_var")),
49 _other_var_num(coupled(
"other_var")),
50 _var_is_pp(getParam<bool>(
"var_is_porepressure")),
51 _viscosity(getParam<Real>(
"fluid_viscosity")),
52 _permeability(getMaterialProperty<RealTensorValue>(
"permeability")),
72 for (
unsigned int nodenum = 0; nodenum <
_num_nodes; ++nodenum)
74 _pp[nodenum] = _var.dofValues()[nodenum];
80 for (
unsigned int nodenum = 0; nodenum <
_num_nodes; ++nodenum)
83 _sat[nodenum] = _var.dofValues()[nodenum];
94 for (
unsigned int nodenum = 0; nodenum <
_num_nodes; ++nodenum)
110 DiracKernel::computeResidual();
116 const Real character =
_character.value(_t, _q_point[_qp]);
117 if (character == 0.0)
120 const Real bh_pressure =
127 const unsigned current_dirac_ptid = currentPointCachedID();
136 if (current_dirac_ptid > 0)
144 _rs[current_dirac_ptid]);
145 if ((character < 0.0 &&
_pp[_i] < bh_pressure) || (character > 0.0 &&
_pp[_i] > bh_pressure))
148 _test[_i][_qp] * std::abs(character) * wc *
_mobility[_i] * (
_pp[_i] - bh_pressure);
151 if (current_dirac_ptid + 1 <
_zs.size() ||
_zs.size() == 1)
158 _rs[current_dirac_ptid]);
159 if ((character < 0.0 &&
_pp[_i] < bh_pressure) || (character > 0.0 &&
_pp[_i] > bh_pressure))
162 _test[_i][_qp] * std::abs(character) * wc *
_mobility[_i] * (
_pp[_i] - bh_pressure);
174 DiracKernel::computeJacobian();
180 return jac(_var.number());
197 const Real character =
_character.value(_t, _q_point[_qp]);
198 if (character == 0.0)
201 const Real bh_pressure =
210 const unsigned current_dirac_ptid = currentPointCachedID();
218 const bool deriv_wrt_pp =
222 if (current_dirac_ptid > 0)
229 _rs[current_dirac_ptid]);
230 if ((character < 0.0 &&
_pp[_i] < bh_pressure) || (character > 0.0 &&
_pp[_i] > bh_pressure))
234 outflowp += std::abs(character) * wc *
237 outflowp += std::abs(character) * wc *
_dmobility_ds[_i] * phi * (
_pp[_i] - bh_pressure);
241 if (current_dirac_ptid <
_zs.size() - 1 ||
_zs.size() == 1)
248 _rs[current_dirac_ptid]);
249 if ((character < 0.0 &&
_pp[_i] < bh_pressure) || (character > 0.0 &&
_pp[_i] > bh_pressure))
253 outflowp += std::abs(character) * wc *
256 outflowp += std::abs(character) * wc *
_dmobility_ds[_i] * phi * (
_pp[_i] - bh_pressure);
260 return _test[_i][_qp] * outflowp;
Approximates a borehole by a sequence of Dirac Points.
InputParameters validParams< Q2PBorehole >()
std::vector< Real > _zs
z points of borehole
virtual Real computeQpResidual()
Computes the Qp residual.
Base class for Richards relative permeability classes that provide relative permeability as a functio...
std::vector< Real > _rs
radii of the borehole
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 > _mobility
nodal mobility
Q2PBorehole(const InputParameters ¶meters)
Creates a new Q2PBorehole This sets all the variables, but also reads the file containing the lines o...
const VariableValue & _other_var_nodal
the other variable in the 2-phase system (this is saturation if Variable=porepressure,...
const unsigned int _other_var_num
the variable number of the other variable
const Real _p_bot
bottomhole pressure of borehole
virtual void computeResidual()
Computes the residual.
std::vector< Real > _dmobility_dp
nodal d(mobility)/d(porepressure)
registerMooseObject("RichardsApp", Q2PBorehole)
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...
virtual void computeJacobian()
Computes the Jacobian.
std::vector< RealTensorValue > _rot_matrix
rotation matrix used in well_constant calculation
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...
Point _bottom_point
the bottom point of the borehole (where bottom_pressure is defined)
InputParameters validParams< PeacemanBorehole >()
const std::string density
Approximates a borehole by a sequence of Dirac Points.
void prepareNodalValues()
calculates the nodal values of pressure, mobility, and derivatives thereof
std::vector< Real > _dmobility_ds
nodal d(mobility)/d(saturation)
const RichardsDensity & _density
fluid density
virtual Real computeQpJacobian()
Computes the diagonal part of the jacobian.
std::vector< Real > _half_seg_len
0.5*(length of polyline segments between points)
Base class for fluid density as a function of porepressure The functions density, ddensity and d2dens...
const Real _viscosity
viscosity
std::vector< Real > _sat
nodal saturation
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 drelperm(Real seff) const =0
derivative of relative permeability wrt effective saturation This must be over-ridden in your derived...
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 bool _var_is_pp
whether the Variable for this BC is porepressure or not
const MaterialProperty< RealTensorValue > & _permeability
permeability
Real jac(unsigned int jvar)
Calculates Jacobian.
const Function & _character
If positive then the borehole acts as a sink (producion well) for porepressure > borehole pressure,...
const RealVectorValue _unit_weight
unit weight of fluid in borehole (for calculating bottomhole pressure at each Dirac Point)
virtual Real computeQpOffDiagJacobian(unsigned int jvar)
Computes the off-diagonal part of the jacobian.
const RichardsRelPerm & _relperm
fluid relative permeability
std::vector< Real > _pp
nodal porepressure
unsigned int _num_nodes
number of nodes in this element.