www.mooseframework.org
RichardsBorehole.h
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 #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 
19 class RichardsBorehole;
20 
21 template <>
22 InputParameters validParams<RichardsBorehole>();
23 
28 {
29 public:
39  RichardsBorehole(const InputParameters & parameters);
40 
46  virtual void computeResidual();
47 
51  virtual Real computeQpResidual();
52 
58  virtual void computeJacobian();
59 
63  virtual Real computeQpJacobian();
64 
71  virtual Real computeQpOffDiagJacobian(unsigned int jvar);
72 
73 protected:
75  const bool _fully_upwind;
76 
79 
81  const unsigned int _num_p;
82 
84  const unsigned int _pvar;
85 
88 
91 
94 
96  unsigned int _num_nodes;
97 
102  std::vector<Real> _mobility;
103 
108  std::vector<std::vector<Real>> _dmobility_dv;
109 
111  const MaterialProperty<std::vector<Real>> & _pp;
112 
114  const MaterialProperty<std::vector<std::vector<Real>>> & _dpp_dv;
115 
117  const MaterialProperty<std::vector<Real>> & _viscosity;
118 
120  const MaterialProperty<RealTensorValue> & _permeability;
121 
123  const MaterialProperty<std::vector<std::vector<Real>>> & _dseff_dv;
124 
126  const MaterialProperty<std::vector<Real>> & _rel_perm;
127 
129  const MaterialProperty<std::vector<std::vector<Real>>> & _drel_perm_dv;
130 
132  const MaterialProperty<std::vector<Real>> & _density;
133 
135  const MaterialProperty<std::vector<std::vector<Real>>> & _ddensity_dv;
136 
145  std::vector<const VariableValue *> _ps_at_nodes;
146 
148  void prepareNodalValues();
149 
154  Real jac(unsigned int wrt_num);
155 };
156 
RichardsBorehole
Approximates a borehole by a sequence of Dirac Points.
Definition: RichardsBorehole.h:27
RichardsBorehole::_ps_at_nodes
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: _...
Definition: RichardsBorehole.h:145
RichardsBorehole::computeResidual
virtual void computeResidual()
Computes the residual.
Definition: RichardsBorehole.C:128
RichardsRelPerm
Base class for Richards relative permeability classes that provide relative permeability as a functio...
Definition: RichardsRelPerm.h:23
RichardsBorehole::_pp
const MaterialProperty< std::vector< Real > > & _pp
fluid porepressure (or porepressures in case of multiphase)
Definition: RichardsBorehole.h:111
RichardsBorehole::_dpp_dv
const MaterialProperty< std::vector< std::vector< Real > > > & _dpp_dv
d(porepressure_i)/d(variable_j)
Definition: RichardsBorehole.h:114
validParams< RichardsBorehole >
InputParameters validParams< RichardsBorehole >()
Definition: RichardsBorehole.C:17
RichardsBorehole::_num_p
const unsigned int _num_p
number of richards variables
Definition: RichardsBorehole.h:81
RichardsVarNames
This holds maps between pressure_var or pressure_var, sat_var used in RichardsMaterial and kernels,...
Definition: RichardsVarNames.h:25
RichardsBorehole::_permeability
const MaterialProperty< RealTensorValue > & _permeability
material permeability
Definition: RichardsBorehole.h:120
RichardsSeff
Base class for effective saturation as a function of porepressure(s) The functions seff,...
Definition: RichardsSeff.h:23
RichardsBorehole::_drel_perm_dv
const MaterialProperty< std::vector< std::vector< Real > > > & _drel_perm_dv
d(relperm_i)/d(variable_j)
Definition: RichardsBorehole.h:129
RichardsBorehole::computeJacobian
virtual void computeJacobian()
Computes the Jacobian.
Definition: RichardsBorehole.C:204
RichardsBorehole::RichardsBorehole
RichardsBorehole(const InputParameters &parameters)
Creates a new RichardsBorehole This sets all the variables, but also reads the file containing the li...
Definition: RichardsBorehole.C:41
RichardsBorehole::_density
const MaterialProperty< std::vector< Real > > & _density
fluid density
Definition: RichardsBorehole.h:132
RichardsBorehole::_density_UO
const RichardsDensity * _density_UO
user object defining the density. Only used if _fully_upwind = true
Definition: RichardsBorehole.h:87
RichardsBorehole::_pvar
const unsigned int _pvar
The moose internal variable number of the richards variable of this Dirac Kernel.
Definition: RichardsBorehole.h:84
RichardsBorehole::_dmobility_dv
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...
Definition: RichardsBorehole.h:108
RichardsVarNames.h
RichardsBorehole::_dseff_dv
const MaterialProperty< std::vector< std::vector< Real > > > & _dseff_dv
deriviatves of Seff wrt variables
Definition: RichardsBorehole.h:123
RichardsBorehole::computeQpJacobian
virtual Real computeQpJacobian()
Computes the diagonal part of the jacobian.
Definition: RichardsBorehole.C:212
RichardsBorehole::_richards_name_UO
const RichardsVarNames & _richards_name_UO
Defines the richards variables in the simulation.
Definition: RichardsBorehole.h:78
PeacemanBorehole
Approximates a borehole by a sequence of Dirac Points.
Definition: PeacemanBorehole.h:25
PeacemanBorehole.h
RichardsBorehole::computeQpResidual
virtual Real computeQpResidual()
Computes the Qp residual.
Definition: RichardsBorehole.C:136
RichardsBorehole::_seff_UO
const RichardsSeff * _seff_UO
user object defining the effective saturation. Only used if _fully_upwind = true
Definition: RichardsBorehole.h:90
RichardsSeff.h
RichardsBorehole::_mobility
std::vector< Real > _mobility
nodal values of mobility = density*relperm/viscosity These are used if _fully_upwind = true
Definition: RichardsBorehole.h:102
RichardsDensity
Base class for fluid density as a function of porepressure The functions density, ddensity and d2dens...
Definition: RichardsDensity.h:24
RichardsBorehole::prepareNodalValues
void prepareNodalValues()
calculates the nodal values of pressure, mobility, and derivatives thereof
Definition: RichardsBorehole.C:88
RichardsBorehole::_ddensity_dv
const MaterialProperty< std::vector< std::vector< Real > > > & _ddensity_dv
d(density_i)/d(variable_j)
Definition: RichardsBorehole.h:135
RichardsRelPerm.h
RichardsBorehole::_fully_upwind
const bool _fully_upwind
Whether to use full upwinding.
Definition: RichardsBorehole.h:75
RichardsBorehole::_relperm_UO
const RichardsRelPerm * _relperm_UO
user object defining the relative permeability. Only used if _fully_upwind = true
Definition: RichardsBorehole.h:93
RichardsBorehole::jac
Real jac(unsigned int wrt_num)
Calculates Jacobian.
Definition: RichardsBorehole.C:230
RichardsDensity.h
RichardsBorehole::_viscosity
const MaterialProperty< std::vector< Real > > & _viscosity
fluid viscosity
Definition: RichardsBorehole.h:117
RichardsBorehole::computeQpOffDiagJacobian
virtual Real computeQpOffDiagJacobian(unsigned int jvar)
Computes the off-diagonal part of the jacobian Note: at March2014 this is never called since moose do...
Definition: RichardsBorehole.C:221
RichardsBorehole::_rel_perm
const MaterialProperty< std::vector< Real > > & _rel_perm
relative permeability
Definition: RichardsBorehole.h:126
RichardsBorehole::_num_nodes
unsigned int _num_nodes
number of nodes in this element. Only used if _fully_upwind = true
Definition: RichardsBorehole.h:96