https://mooseframework.inl.gov
RichardsBorehole.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
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 #pragma once
11 
12 // Moose Includes
13 #include "PeacemanBorehole.h"
14 #include "RichardsVarNames.h"
15 #include "RichardsDensity.h"
16 #include "RichardsRelPerm.h"
17 #include "RichardsSeff.h"
18 
23 {
24 public:
35 
37 
43  virtual void computeResidual();
44 
48  virtual Real computeQpResidual();
49 
55  virtual void computeJacobian();
56 
60  virtual Real computeQpJacobian();
61 
68  virtual Real computeQpOffDiagJacobian(unsigned int jvar);
69 
70 protected:
72  const bool _fully_upwind;
73 
76 
78  const unsigned int _num_p;
79 
81  const unsigned int _pvar;
82 
85 
87  const RichardsSeff * const _seff_UO;
88 
91 
93  unsigned int _num_nodes;
94 
99  std::vector<Real> _mobility;
100 
105  std::vector<std::vector<Real>> _dmobility_dv;
106 
109 
112 
115 
118 
121 
124 
127 
130 
133 
142  std::vector<const VariableValue *> _ps_at_nodes;
143 
145  void prepareNodalValues();
146 
151  Real jac(unsigned int wrt_num);
152 };
const MaterialProperty< RealTensorValue > & _permeability
material permeability
const bool _fully_upwind
Whether to use full upwinding.
const MaterialProperty< std::vector< Real > > & _density
fluid density
Real jac(unsigned int wrt_num)
Calculates Jacobian.
const MaterialProperty< std::vector< Real > > & _viscosity
fluid viscosity
Base class for effective saturation as a function of porepressure(s) The functions seff...
Definition: RichardsSeff.h:18
Base class for Richards relative permeability classes that provide relative permeability as a functio...
const MaterialProperty< std::vector< std::vector< Real > > > & _dseff_dv
deriviatves of Seff wrt variables
virtual Real computeQpJacobian()
Computes the diagonal part of the jacobian.
unsigned int _num_nodes
number of nodes in this element. Only used if _fully_upwind = true
Approximates a borehole by a sequence of Dirac Points.
This holds maps between pressure_var or pressure_var, sat_var used in RichardsMaterial and kernels...
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 RichardsRelPerm *const _relperm_UO
user object defining the relative permeability. Only used if _fully_upwind = true ...
virtual void computeResidual()
Computes the residual.
const RichardsDensity *const _density_UO
user object defining the density. Only used if _fully_upwind = true
static InputParameters validParams()
Creates a new RichardsBorehole This sets all the variables, but also reads the file containing the li...
const MaterialProperty< std::vector< Real > > & _rel_perm
relative permeability
const RichardsSeff *const _seff_UO
user object defining the effective saturation. Only used if _fully_upwind = true
virtual Real computeQpOffDiagJacobian(unsigned int jvar)
Computes the off-diagonal part of the jacobian Note: at March2014 this is never called since moose do...
const MaterialProperty< std::vector< std::vector< Real > > > & _drel_perm_dv
d(relperm_i)/d(variable_j)
virtual void computeJacobian()
Computes the Jacobian.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
RichardsBorehole(const InputParameters &parameters)
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 InputParameters & parameters() const
Base class for fluid density as a function of porepressure The functions density, ddensity and d2dens...
const RichardsVarNames & _richards_name_UO
Defines the richards variables in the simulation.
virtual Real computeQpResidual()
Computes the Qp residual.
void prepareNodalValues()
calculates the nodal values of pressure, mobility, and derivatives thereof
const unsigned int _num_p
number of richards variables
const MaterialProperty< std::vector< std::vector< Real > > > & _ddensity_dv
d(density_i)/d(variable_j)
Approximates a borehole by a sequence of Dirac Points.