www.mooseframework.org
RichardsMaterial.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 #include "Material.h"
13 
14 #include "RichardsVarNames.h"
15 #include "RichardsDensity.h"
16 #include "RichardsRelPerm.h"
17 #include "RichardsSeff.h"
18 #include "RichardsSat.h"
19 #include "RichardsSUPG.h"
20 
21 // Forward Declarations
22 
23 class RichardsMaterial : public Material
24 {
25 public:
27 
29 
30 protected:
33 
37 
40 
43 
49 
52  unsigned int _num_p;
53 
54  std::vector<const RichardsRelPerm *> _material_relperm_UO;
55  std::vector<const RichardsSeff *> _material_seff_UO;
56  std::vector<const RichardsSat *> _material_sat_UO;
57  std::vector<const RichardsDensity *> _material_density_UO;
58  std::vector<const RichardsSUPG *> _material_SUPG_UO;
59 
60  const std::vector<const VariableValue *> _perm_change;
61 
62  virtual void computeProperties();
63 
71  void computePandSeff();
72 
79  void computeDerivedQuantities(unsigned int qp);
80 
86  void compute2ndDerivedQuantities(unsigned int qp);
87 
93  void zeroSUPG(unsigned int qp);
94 
96  void computeSUPG();
97 
98 private:
101 
102  std::vector<Real> _material_viscosity;
103 
106 
109 
112 
115 
118 
121 
124 
127 
130 
132  MaterialProperty<std::vector<Real>> & _seff; // effective saturation
133 
136 
139 
142 
145 
148 
151 
154 
157 
160 
163 
166 
169 
172 
175 
178 
181 
184 
187 
190 
193  _dtauvel_SUPG_dgradp; // d (_tauvel_SUPG_i)/d(_grad_variable_j)
195  _dtauvel_SUPG_dp; // d (_tauvel_SUPG_i)/d(variable_j)
196 
198  std::vector<std::vector<std::vector<Real>>> _d2density;
199 
201  std::vector<std::vector<std::vector<Real>>> _d2rel_perm_dv;
202 
203  std::vector<const VariableValue *> _pressure_vals;
204  std::vector<const VariableValue *> _pressure_old_vals;
205  std::vector<const VariableGradient *> _grad_p;
206 
214  void zero2ndDerivedQuantities(unsigned int qp);
215 };
MaterialProperty< std::vector< std::vector< RealTensorValue > > > & _dflux_no_mob_dgradv
d(_flux_no_mob_i)/d(grad(variable_j))
MaterialProperty< std::vector< std::vector< std::vector< Real > > > > & _d2pp_dv
d^2(porepressure_i)/d(variable_j)/d(variable_k)
MaterialProperty< std::vector< Real > > & _pp
porepressure(s)
MaterialProperty< RealVectorValue > & _gravity
MaterialProperty< std::vector< std::vector< std::vector< RealTensorValue > > > > & _d2flux_dvdgradv
d^2(Richards flux_i)/d(variable_j)/d(grad(variable_k)), here flux_i is the i_th flux, which is itself a RealVectorValue. We should have _d2flux_dvdgradv[i][j][k] = _d2flux_dgradvdv[i][k][j], but i think it is more clear having both, and hopefully not a blowout on memory/CPU.
std::vector< std::vector< std::vector< Real > > > _d2rel_perm_dv
d^2(relperm_i)/dP_j/dP_k - used in various derivative calculations
RichardsMaterial(const InputParameters &parameters)
MaterialProperty< std::vector< std::vector< std::vector< RealTensorValue > > > > & _d2flux_dgradvdv
d^2(Richards flux_i)/d(grad(variable_j))/d(variable_k), here flux_i is the i_th flux, which is itself a RealVectorValue
std::vector< const RichardsRelPerm * > _material_relperm_UO
MaterialProperty< std::vector< Real > > & _pp_old
old values of porepressure(s)
const RichardsVarNames & _richards_name_UO
The variable names userobject for the Richards variables.
MaterialProperty< std::vector< RealVectorValue > > & _flux
fluid flux (a vector of fluxes for multicomponent)
std::vector< Real > _material_viscosity
const std::vector< const VariableValue * > _perm_change
unsigned int _num_p
This holds maps between pressure_var or pressure_var, sat_var used in RichardsMaterial and kernels...
MaterialProperty< std::vector< Real > > & _sat
saturation (vector of saturations in case of multiphase)
std::vector< const RichardsSUPG * > _material_SUPG_UO
MaterialProperty< std::vector< std::vector< Real > > > & _drel_perm_dv
d(relperm_i)/d(variable_j)
Real _material_por
porosity as entered by the user
const VariableValue & _por_change
porosity changes. if not entered they default to zero
void zero2ndDerivedQuantities(unsigned int qp)
Zeroes 2nd derivatives of the flux.
TensorValue< Real > RealTensorValue
MaterialProperty< std::vector< Real > > & _seff_old
old effective saturation
MaterialProperty< std::vector< std::vector< Real > > > & _dpp_dv
d(porepressure_i)/d(variable_j)
void computePandSeff()
computes the quadpoint values of the porepressure(s) and effective saturation(s), and their derivativ...
MaterialProperty< std::vector< std::vector< std::vector< RealVectorValue > > > > & _d2flux_dvdv
d^2(Richards flux_i)/d(variable_j)/d(variable_k), here flux_i is the i_th flux, which is itself a Rea...
MaterialProperty< RealTensorValue > & _permeability
MaterialProperty< std::vector< std::vector< RealVectorValue > > > & _dflux_no_mob_dv
d(_flux_no_mob_i)/d(variable_j)
MaterialProperty< std::vector< std::vector< RealTensorValue > > > & _dtauvel_SUPG_dgradp
std::vector< const VariableGradient * > _grad_p
MaterialProperty< std::vector< std::vector< RealVectorValue > > > & _dflux_dv
d(Richards flux_i)/d(variable_j), here flux_i is the i_th flux, which is itself a RealVectorValue ...
std::vector< const VariableValue * > _pressure_vals
MaterialProperty< std::vector< Real > > & _rel_perm
relative permeability (vector of relative permeabilities in case of multiphase)
MaterialProperty< std::vector< Real > > & _mass
fluid mass (a vector of masses for multicomponent)
RealVectorValue _material_gravity
gravity as entered by user
MaterialProperty< std::vector< RealVectorValue > > & _flux_no_mob
permeability*(grad(P) - density*gravity) (a vector of these for multicomponent)
std::vector< const RichardsDensity * > _material_density_UO
void computeDerivedQuantities(unsigned int qp)
Computes the "derived" quantities — those that depend on porepressure(s) and effective saturation(s)...
std::vector< const VariableValue * > _pressure_old_vals
MaterialProperty< std::vector< Real > > & _density
fluid density (or densities for multiphase problems)
MaterialProperty< std::vector< Real > > & _sat_old
old saturation
OutputTools< Real >::VariableValue VariableValue
std::vector< const RichardsSeff * > _material_seff_UO
MaterialProperty< std::vector< Real > > & _mass_old
old value of fluid mass (a vector of masses for multicomponent)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
MaterialProperty< Real > & _porosity
void computeSUPG()
Computes the tauvel_SUPG and its derivatives.
MaterialProperty< std::vector< RealVectorValue > > & _tauvel_SUPG
MaterialProperty< std::vector< std::vector< std::vector< Real > > > > & _d2seff_dv
d^2(Seff_i)/d(variable_j)/d(variable_k)
void compute2ndDerivedQuantities(unsigned int qp)
Computes 2nd derivatives of the flux.
std::vector< const RichardsSat * > _material_sat_UO
void zeroSUPG(unsigned int qp)
Assigns and zeroes the MaterialProperties associated with SUPG.
MaterialProperty< std::vector< Real > > & _density_old
old fluid density (or densities for multiphase problems)
MaterialProperty< Real > & _porosity_old
material properties
const VariableValue & _por_change_old
MaterialProperty< std::vector< Real > > & _seff
effective saturation (vector of effective saturations in case of multiphase)
RealTensorValue _material_perm
permeability as entered by the user
const InputParameters & parameters() const
MaterialProperty< std::vector< std::vector< RealTensorValue > > > & _dflux_dgradv
d(Richards flux_i)/d(grad(variable_j)), here flux_i is the i_th flux, which is itself a RealVectorVal...
Real _trace_perm
trace of permeability tensor
MaterialProperty< std::vector< std::vector< Real > > > & _ddensity_dv
d(density_i)/d(variable_j)
MaterialProperty< std::vector< std::vector< Real > > > & _dsat_dv
d(saturation_i)/d(variable_j)
MaterialProperty< std::vector< std::vector< Real > > > & _dseff_dv
d(Seff_i)/d(variable_j)
MaterialProperty< std::vector< std::vector< Real > > > & _dmass
d(fluid mass_i)/dP_j (a vector of masses for multicomponent)
static InputParameters validParams()
MaterialProperty< std::vector< std::vector< RealVectorValue > > > & _dtauvel_SUPG_dp
std::vector< std::vector< std::vector< Real > > > _d2density
d^2(density)/dp_j/dP_k - used in various derivative calculations
MaterialProperty< std::vector< Real > > & _viscosity
fluid viscosity (or viscosities in the multiphase case)
virtual void computeProperties()