www.mooseframework.org
Q2PBorehole.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #include "Q2PBorehole.h"
11 #include "RotationMatrix.h"
12 
13 registerMooseObject("RichardsApp", Q2PBorehole);
14 
15 template <>
16 InputParameters
18 {
19  InputParameters params = validParams<PeacemanBorehole>();
20  params.addRequiredParam<UserObjectName>(
21  "fluid_density",
22  "A RichardsDensity UserObject that defines the fluid density as a function of pressure.");
23  params.addRequiredParam<UserObjectName>(
24  "fluid_relperm",
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 "
30  "vice-versa.");
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");
41  return params;
42 }
43 
44 Q2PBorehole::Q2PBorehole(const InputParameters & parameters)
45  : PeacemanBorehole(parameters),
46  _density(getUserObject<RichardsDensity>("fluid_density")),
47  _relperm(getUserObject<RichardsRelPerm>("fluid_relperm")),
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")),
53  _num_nodes(0),
54  _pp(0),
55  _sat(0),
56  _mobility(0),
57  _dmobility_dp(0),
58  _dmobility_ds(0)
59 {
60 }
61 
62 void
64 {
66 
67  // set _pp and _sat variables
68  _pp.resize(_num_nodes);
69  _sat.resize(_num_nodes);
70  if (_var_is_pp)
71  {
72  for (unsigned int nodenum = 0; nodenum < _num_nodes; ++nodenum)
73  {
74  _pp[nodenum] = _var.dofValues()[nodenum];
75  _sat[nodenum] = _other_var_nodal[nodenum];
76  }
77  }
78  else
79  {
80  for (unsigned int nodenum = 0; nodenum < _num_nodes; ++nodenum)
81  {
82  _pp[nodenum] = _other_var_nodal[nodenum];
83  _sat[nodenum] = _var.dofValues()[nodenum];
84  }
85  }
86 
87  Real density;
88  Real ddensity_dp;
89  Real relperm;
90  Real drelperm_ds;
91  _mobility.resize(_num_nodes);
92  _dmobility_dp.resize(_num_nodes);
93  _dmobility_ds.resize(_num_nodes);
94  for (unsigned int nodenum = 0; nodenum < _num_nodes; ++nodenum)
95  {
96  density = _density.density(_pp[nodenum]);
97  ddensity_dp = _density.ddensity(_pp[nodenum]);
98  relperm = _relperm.relperm(_sat[nodenum]);
99  drelperm_ds = _relperm.drelperm(_sat[nodenum]);
100  _mobility[nodenum] = density * relperm / _viscosity;
101  _dmobility_dp[nodenum] = ddensity_dp * relperm / _viscosity;
102  _dmobility_ds[nodenum] = density * drelperm_ds / _viscosity;
103  }
104 }
105 
106 void
108 {
110  DiracKernel::computeResidual();
111 }
112 
113 Real
115 {
116  const Real character = _character.value(_t, _q_point[_qp]);
117  if (character == 0.0)
118  return 0.0;
119 
120  const Real bh_pressure =
121  _p_bot +
122  _unit_weight *
123  (_q_point[_qp] -
124  _bottom_point); // really want to use _q_point instaed of _current_point, i think?!
125 
126  // Get the ID we initially assigned to this point
127  const unsigned current_dirac_ptid = currentPointCachedID();
128 
129  // If getting the ID failed, fall back to the old bodge!
130  // if (current_dirac_ptid == libMesh::invalid_uint)
131  // current_dirac_ptid = (_zs.size() > 2) ? 1 : 0;
132 
133  Real outflow(0.0); // this is the flow rate from porespace out of the system
134 
135  Real wc(0.0);
136  if (current_dirac_ptid > 0)
137  // contribution from half-segment "behind" this point (must have >1 point for
138  // current_dirac_ptid>0)
139  {
140  wc = wellConstant(_permeability[0],
141  _rot_matrix[current_dirac_ptid - 1],
142  _half_seg_len[current_dirac_ptid - 1],
143  _current_elem,
144  _rs[current_dirac_ptid]);
145  if ((character < 0.0 && _pp[_i] < bh_pressure) || (character > 0.0 && _pp[_i] > bh_pressure))
146  // injection, so outflow<0 || // production, so outflow>0
147  outflow +=
148  _test[_i][_qp] * std::abs(character) * wc * _mobility[_i] * (_pp[_i] - bh_pressure);
149  }
150 
151  if (current_dirac_ptid + 1 < _zs.size() || _zs.size() == 1)
152  // contribution from half-segment "ahead of" this point, or we only have one point
153  {
154  wc = wellConstant(_permeability[0],
155  _rot_matrix[current_dirac_ptid],
156  _half_seg_len[current_dirac_ptid],
157  _current_elem,
158  _rs[current_dirac_ptid]);
159  if ((character < 0.0 && _pp[_i] < bh_pressure) || (character > 0.0 && _pp[_i] > bh_pressure))
160  // injection, so outflow<0 || // production, so outflow>0
161  outflow +=
162  _test[_i][_qp] * std::abs(character) * wc * _mobility[_i] * (_pp[_i] - bh_pressure);
163  }
164 
166  outflow * _dt); // this is not thread safe, but DiracKernel's aren't currently threaded
167  return outflow;
168 }
169 
170 void
172 {
174  DiracKernel::computeJacobian();
175 }
176 
177 Real
179 {
180  return jac(_var.number());
181 }
182 
183 Real
185 {
186  if (jvar == _other_var_num || jvar == _var.number())
187  return jac(jvar);
188  return 0.0;
189 }
190 
191 Real
192 Q2PBorehole::jac(unsigned int jvar)
193 {
194  if (_i != _j)
195  return 0.0;
196 
197  const Real character = _character.value(_t, _q_point[_qp]);
198  if (character == 0.0)
199  return 0.0;
200 
201  const Real bh_pressure =
202  _p_bot +
203  _unit_weight *
204  (_q_point[_qp] -
205  _bottom_point); // really want to use _q_point instaed of _current_point, i think?!
206 
207  const Real phi = 1;
208 
209  // Get the ID we initially assigned to this point
210  const unsigned current_dirac_ptid = currentPointCachedID();
211 
212  // If getting the ID failed, fall back to the old bodge!
213  // if (current_dirac_ptid == libMesh::invalid_uint)
214  // current_dirac_ptid = (_zs.size() > 2) ? 1 : 0;
215 
216  Real outflowp(0.0);
217 
218  const bool deriv_wrt_pp =
219  (_var_is_pp && (jvar == _var.number())) || (!_var_is_pp && (jvar == _other_var_num));
220 
221  Real wc(0.0);
222  if (current_dirac_ptid > 0)
223  // contribution from half-segment "behind" this point
224  {
225  wc = wellConstant(_permeability[0],
226  _rot_matrix[current_dirac_ptid - 1],
227  _half_seg_len[current_dirac_ptid - 1],
228  _current_elem,
229  _rs[current_dirac_ptid]);
230  if ((character < 0.0 && _pp[_i] < bh_pressure) || (character > 0.0 && _pp[_i] > bh_pressure))
231  {
232  // injection, so outflow<0 || // production, so outflow>0
233  if (deriv_wrt_pp)
234  outflowp += std::abs(character) * wc *
235  (_mobility[_i] * phi + _dmobility_dp[_i] * phi * (_pp[_i] - bh_pressure));
236  else
237  outflowp += std::abs(character) * wc * _dmobility_ds[_i] * phi * (_pp[_i] - bh_pressure);
238  }
239  }
240 
241  if (current_dirac_ptid < _zs.size() - 1 || _zs.size() == 1)
242  // contribution from half-segment "ahead of" this point
243  {
244  wc = wellConstant(_permeability[0],
245  _rot_matrix[current_dirac_ptid],
246  _half_seg_len[current_dirac_ptid],
247  _current_elem,
248  _rs[current_dirac_ptid]);
249  if ((character < 0.0 && _pp[_i] < bh_pressure) || (character > 0.0 && _pp[_i] > bh_pressure))
250  {
251  // injection, so outflow<0 || // production, so outflow>0
252  if (deriv_wrt_pp)
253  outflowp += std::abs(character) * wc *
254  (_mobility[_i] * phi + _dmobility_dp[_i] * phi * (_pp[_i] - bh_pressure));
255  else
256  outflowp += std::abs(character) * wc * _dmobility_ds[_i] * phi * (_pp[_i] - bh_pressure);
257  }
258  }
259 
260  return _test[_i][_qp] * outflowp;
261 }
Q2PBorehole
Approximates a borehole by a sequence of Dirac Points.
Definition: Q2PBorehole.h:26
validParams< Q2PBorehole >
InputParameters validParams< Q2PBorehole >()
Definition: Q2PBorehole.C:17
PeacemanBorehole::_zs
std::vector< Real > _zs
z points of borehole
Definition: PeacemanBorehole.h:91
Q2PBorehole::computeQpResidual
virtual Real computeQpResidual()
Computes the Qp residual.
Definition: Q2PBorehole.C:114
RichardsRelPerm
Base class for Richards relative permeability classes that provide relative permeability as a functio...
Definition: RichardsRelPerm.h:23
PeacemanBorehole::_rs
std::vector< Real > _rs
radii of the borehole
Definition: PeacemanBorehole.h:82
RichardsRelPerm::relperm
virtual Real relperm(Real seff) const =0
relative permeability as a function of effective saturation This must be over-ridden in your derived ...
Q2PBorehole::_mobility
std::vector< Real > _mobility
nodal mobility
Definition: Q2PBorehole.h:101
Q2PBorehole::Q2PBorehole
Q2PBorehole(const InputParameters &parameters)
Creates a new Q2PBorehole This sets all the variables, but also reads the file containing the lines o...
Definition: Q2PBorehole.C:44
Q2PBorehole::_other_var_nodal
const VariableValue & _other_var_nodal
the other variable in the 2-phase system (this is saturation if Variable=porepressure,...
Definition: Q2PBorehole.h:77
Q2PBorehole::_other_var_num
const unsigned int _other_var_num
the variable number of the other variable
Definition: Q2PBorehole.h:80
PeacemanBorehole::_p_bot
const Real _p_bot
bottomhole pressure of borehole
Definition: PeacemanBorehole.h:69
Q2PBorehole::computeResidual
virtual void computeResidual()
Computes the residual.
Definition: Q2PBorehole.C:107
Q2PBorehole::_dmobility_dp
std::vector< Real > _dmobility_dp
nodal d(mobility)/d(porepressure)
Definition: Q2PBorehole.h:104
registerMooseObject
registerMooseObject("RichardsApp", Q2PBorehole)
PeacemanBorehole::wellConstant
Real wellConstant(const RealTensorValue &perm, const RealTensorValue &rot, const Real &half_len, const Elem *ele, const Real &rad)
Calculates Peaceman's form of the borehole well constant Z Chen, Y Zhang, Well flow models for variou...
Definition: PeacemanBorehole.C:194
Q2PBorehole.h
Q2PBorehole::computeJacobian
virtual void computeJacobian()
Computes the Jacobian.
Definition: Q2PBorehole.C:171
PeacemanBorehole::_rot_matrix
std::vector< RealTensorValue > _rot_matrix
rotation matrix used in well_constant calculation
Definition: PeacemanBorehole.h:100
RichardsDensity::density
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...
PeacemanBorehole::_bottom_point
Point _bottom_point
the bottom point of the borehole (where bottom_pressure is defined)
Definition: PeacemanBorehole.h:94
validParams< PeacemanBorehole >
InputParameters validParams< PeacemanBorehole >()
Definition: PeacemanBorehole.C:17
NS::density
const std::string density
Definition: NS.h:16
PeacemanBorehole
Approximates a borehole by a sequence of Dirac Points.
Definition: PeacemanBorehole.h:25
Q2PBorehole::prepareNodalValues
void prepareNodalValues()
calculates the nodal values of pressure, mobility, and derivatives thereof
Definition: Q2PBorehole.C:63
Q2PBorehole::_dmobility_ds
std::vector< Real > _dmobility_ds
nodal d(mobility)/d(saturation)
Definition: Q2PBorehole.h:107
Q2PBorehole::_density
const RichardsDensity & _density
fluid density
Definition: Q2PBorehole.h:71
Q2PBorehole::computeQpJacobian
virtual Real computeQpJacobian()
Computes the diagonal part of the jacobian.
Definition: Q2PBorehole.C:178
PeacemanBorehole::_half_seg_len
std::vector< Real > _half_seg_len
0.5*(length of polyline segments between points)
Definition: PeacemanBorehole.h:97
RichardsDensity
Base class for fluid density as a function of porepressure The functions density, ddensity and d2dens...
Definition: RichardsDensity.h:24
Q2PBorehole::_viscosity
const Real _viscosity
viscosity
Definition: Q2PBorehole.h:86
Q2PBorehole::_sat
std::vector< Real > _sat
nodal saturation
Definition: Q2PBorehole.h:98
RichardsDensity::ddensity
virtual Real ddensity(Real p) const =0
derivative of fluid density wrt porepressure This must be over-ridden in derived classes to provide a...
RichardsRelPerm::drelperm
virtual Real drelperm(Real seff) const =0
derivative of relative permeability wrt effective saturation This must be over-ridden in your derived...
RichardsSumQuantity::add
void add(Real contrib)
adds contrib to _total
Definition: RichardsSumQuantity.C:37
PeacemanBorehole::_total_outflow_mass
RichardsSumQuantity & _total_outflow_mass
This is used to hold the total fluid flowing into the borehole Hence, it is positive for production w...
Definition: PeacemanBorehole.h:79
Q2PBorehole::_var_is_pp
const bool _var_is_pp
whether the Variable for this BC is porepressure or not
Definition: Q2PBorehole.h:83
Q2PBorehole::_permeability
const MaterialProperty< RealTensorValue > & _permeability
permeability
Definition: Q2PBorehole.h:89
Q2PBorehole::jac
Real jac(unsigned int jvar)
Calculates Jacobian.
Definition: Q2PBorehole.C:192
PeacemanBorehole::_character
const Function & _character
If positive then the borehole acts as a sink (producion well) for porepressure > borehole pressure,...
Definition: PeacemanBorehole.h:66
PeacemanBorehole::_unit_weight
const RealVectorValue _unit_weight
unit weight of fluid in borehole (for calculating bottomhole pressure at each Dirac Point)
Definition: PeacemanBorehole.h:72
Q2PBorehole::computeQpOffDiagJacobian
virtual Real computeQpOffDiagJacobian(unsigned int jvar)
Computes the off-diagonal part of the jacobian.
Definition: Q2PBorehole.C:184
Q2PBorehole::_relperm
const RichardsRelPerm & _relperm
fluid relative permeability
Definition: Q2PBorehole.h:74
Q2PBorehole::_pp
std::vector< Real > _pp
nodal porepressure
Definition: Q2PBorehole.h:95
Q2PBorehole::_num_nodes
unsigned int _num_nodes
number of nodes in this element.
Definition: Q2PBorehole.h:92