https://mooseframework.inl.gov
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
PorousFlowFluidStateTempl< is_ad > Class Template Reference

Fluid state class using a persistent set of primary variables for the mutliphase, multicomponent case. More...

#include <PorousFlowFluidState.h>

Inheritance diagram for PorousFlowFluidStateTempl< is_ad >:
[legend]

Public Types

typedef DerivativeMaterialPropertyNameInterface::SymbolName SymbolName
 

Public Member Functions

 PorousFlowFluidStateTempl (const InputParameters &parameters)
 
const GenericMaterialProperty< U, is_ad > & getDefaultMaterialProperty (const std::string &name)
 
const GenericMaterialProperty< U, is_ad > & getDefaultMaterialPropertyByName (const std::string &name)
 
void validateDerivativeMaterialPropertyBase (const std::string &base)
 
const MaterialPropertyName derivativePropertyName (const MaterialPropertyName &base, const std::vector< SymbolName > &c) const
 
const MaterialPropertyName derivativePropertyNameFirst (const MaterialPropertyName &base, const SymbolName &c1) const
 
const MaterialPropertyName derivativePropertyNameSecond (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2) const
 
const MaterialPropertyName derivativePropertyNameThird (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2, const SymbolName &c3) const
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const std::vector< VariableName > &c)
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const std::vector< SymbolName > &c)
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const std::vector< VariableName > &c)
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const std::vector< SymbolName > &c)
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const std::vector< VariableName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const std::vector< SymbolName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const SymbolName &c1, unsigned int v2, unsigned int v3=libMesh::invalid_uint)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, unsigned int v1, unsigned int v2=libMesh::invalid_uint, unsigned int v3=libMesh::invalid_uint)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const std::vector< VariableName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const std::vector< SymbolName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const SymbolName &c1, unsigned int v2, unsigned int v3=libMesh::invalid_uint)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, unsigned int v1, unsigned int v2=libMesh::invalid_uint, unsigned int v3=libMesh::invalid_uint)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const std::vector< VariableName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const std::vector< SymbolName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const std::vector< VariableName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const std::vector< SymbolName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
void validateCoupling (const MaterialPropertyName &base, const std::vector< VariableName > &c, bool validate_aux=true)
 
void validateCoupling (const MaterialPropertyName &base, const VariableName &c1="", const VariableName &c2="", const VariableName &c3="")
 
void validateCoupling (const MaterialPropertyName &base, const std::vector< VariableName > &c, bool validate_aux=true)
 
void validateCoupling (const MaterialPropertyName &base, const VariableName &c1="", const VariableName &c2="", const VariableName &c3="")
 
void validateNonlinearCoupling (const MaterialPropertyName &base, const VariableName &c1="", const VariableName &c2="", const VariableName &c3="")
 
void validateNonlinearCoupling (const MaterialPropertyName &base, const VariableName &c1="", const VariableName &c2="", const VariableName &c3="")
 
const MaterialPropertyName propertyName (const MaterialPropertyName &base, const std::vector< SymbolName > &c) const
 
const MaterialPropertyName propertyName (const MaterialPropertyName &base, const std::vector< SymbolName > &c) const
 
const MaterialPropertyName propertyNameFirst (const MaterialPropertyName &base, const SymbolName &c1) const
 
const MaterialPropertyName propertyNameFirst (const MaterialPropertyName &base, const SymbolName &c1) const
 
const MaterialPropertyName propertyNameSecond (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2) const
 
const MaterialPropertyName propertyNameSecond (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2) const
 
const MaterialPropertyName propertyNameThird (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2, const SymbolName &c3) const
 
const MaterialPropertyName propertyNameThird (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2, const SymbolName &c3) const
 

Static Public Member Functions

static InputParameters validParams ()
 

Protected Member Functions

virtual void initQpStatefulProperties () override
 
virtual void computeQpProperties () override
 
virtual void thermophysicalProperties () override
 Calculates all required thermophysical properties and derivatives for each phase and fluid component. More...
 
virtual void setMaterialVectorSize () const
 Size material property vectors and initialise with zeros. More...
 
GenericReal< is_ad > genericValue (const ADReal &value)
 Templated method that takes an ADReal value, and returns the raw value when is_ad = false, and the AD value when is_ad = true. More...
 
template<>
GenericReal< false > genericValue (const ADReal &value)
 
template<>
GenericReal< true > genericValue (const ADReal &value)
 

Protected Attributes

const GenericVariableValue< is_ad > & _gas_porepressure
 Porepressure. More...
 
const GenericVariableGradient< is_ad > & _gas_gradp_qp
 Gradient of porepressure (only defined at the qps) More...
 
const unsigned int _gas_porepressure_varnum
 Moose variable number of the gas porepressure. More...
 
const unsigned int _pvar
 PorousFlow variable number of the gas porepressure. More...
 
std::vector< const GenericVariableValue< is_ad > * > _Z
 Total mass fraction(s) of the gas component(s) summed over all phases. More...
 
std::vector< const GenericVariableGradient< is_ad > * > _gradZ_qp
 Gradient(s) of total mass fraction(s) of the gas component(s) (only defined at the qps) More...
 
std::vector< unsigned int_Z_varnum
 Moose variable number of Z. More...
 
std::vector< unsigned int_Zvar
 PorousFlow variable number of Z. More...
 
const unsigned int _num_Z_vars
 Number of coupled total mass fractions. Should be _num_phases - 1. More...
 
const bool _is_Xnacl_nodal
 Flag for nodal NaCl mass fraction. More...
 
const GenericVariableValue< is_ad > & _Xnacl
 Salt mass fraction (kg/kg) More...
 
const GenericVariableGradient< is_ad > & _grad_Xnacl_qp
 Gradient of salt mass fraction (only defined at the qps) More...
 
const unsigned int _Xnacl_varnum
 Salt mass fraction variable number. More...
 
const unsigned int _Xvar
 Salt mass fraction PorousFlow variable number. More...
 
const PorousFlowFluidStateMultiComponentBase_fs
 FluidState UserObject. More...
 
const unsigned int _aqueous_phase_number
 Phase number of the aqueous phase. More...
 
const unsigned int _gas_phase_number
 Phase number of the gas phase. More...
 
const unsigned int _aqueous_fluid_component
 Fluid component number of the aqueous component. More...
 
const unsigned int _gas_fluid_component
 Fluid component number of the gas phase. More...
 
const unsigned int _salt_component
 Salt component index. More...
 
const GenericMaterialProperty< Real, is_ad > & _temperature
 Temperature. More...
 
const GenericMaterialProperty< RealGradient, is_ad > *const _gradT_qp
 Gradient of temperature (only defined at the qps) More...
 
const MaterialProperty< std::vector< Real > > *const _dtemperature_dvar
 Derivative of temperature wrt PorousFlow variables. More...
 
const unsigned int _temperature_varnum
 Moose variable number of the temperature. More...
 
const unsigned int _Tvar
 PorousFlow variable number of the temperature. More...
 
const unsigned int _pidx
 Index of derivative wrt pressure. More...
 
const unsigned int _Tidx
 Index of derivative wrt temperature. More...
 
const unsigned int _Zidx
 Index of derivative wrt total mass fraction Z. More...
 
const unsigned int _Xidx
 Index of derivative wrt salt mass fraction X. More...
 
 usingPorousFlowFluidStateBaseMaterialMembers
 
 usingPorousFlowFluidStateBaseMaterialDerivativeMembers
 
const Real _T_c2k
 Conversion from degrees Celsius to degrees Kelvin. More...
 
bool _is_initqp
 Flag to indicate whether stateful properties should be computed. More...
 
std::vector< FluidStateProperties_fsp
 FluidStateProperties data for each fluid component in each fluid phase. More...
 
const PorousFlowCapillaryPressure_pc
 Capillary pressure UserObject. More...
 
const std::string _sfx
 Suffix to append to material property names (either qp or nodal as required) More...
 
GenericMaterialProperty< std::vector< std::vector< Real > >, is_ad > & _mass_frac
 Mass fraction matrix (indexing is fluid component in fluid phase) More...
 
GenericMaterialProperty< std::vector< std::vector< RealGradient > >, is_ad > *const _grad_mass_frac_qp
 Gradient of the mass fraction matrix (only defined at the qps) More...
 
MaterialProperty< std::vector< std::vector< std::vector< Real > > > > *const _dmass_frac_dvar
 Derivative of the mass fraction matrix with respect to the Porous Flow variables. More...
 
GenericMaterialProperty< std::vector< Real >, is_ad > & _fluid_density
 Fluid density of each phase. More...
 
MaterialProperty< std::vector< std::vector< Real > > > *const _dfluid_density_dvar
 Derivative of the fluid density for each phase wrt PorousFlow variables. More...
 
GenericMaterialProperty< std::vector< Real >, is_ad > & _fluid_viscosity
 Viscosity of each phase. More...
 
MaterialProperty< std::vector< std::vector< Real > > > *const _dfluid_viscosity_dvar
 Derivative of the fluid viscosity for each phase wrt PorousFlow variables. More...
 
GenericMaterialProperty< std::vector< Real >, is_ad > & _fluid_enthalpy
 Enthalpy of each phase. More...
 
MaterialProperty< std::vector< std::vector< Real > > > *const _dfluid_enthalpy_dvar
 Derivative of the fluid enthalpy for each phase wrt PorousFlow variables. More...
 
GenericMaterialProperty< std::vector< Real >, is_ad > & _fluid_internal_energy
 Internal energy of each phase. More...
 
MaterialProperty< std::vector< std::vector< Real > > > *const _dfluid_internal_energy_dvar
 Derivative of the fluid internal energy for each phase wrt PorousFlow variables. More...
 
 usingPorousFlowVariableBaseMembers
 
const unsigned int _num_phases
 Number of phases. More...
 
const unsigned int _num_components
 Number of components. More...
 
const unsigned int _num_pf_vars
 Number of PorousFlow variables. More...
 
GenericMaterialProperty< std::vector< Real >, is_ad > & _porepressure
 Computed nodal or quadpoint values of porepressure of the phases. More...
 
MaterialProperty< std::vector< std::vector< Real > > > *const _dporepressure_dvar
 d(porepressure)/d(PorousFlow variable) More...
 
GenericMaterialProperty< std::vector< RealGradient >, is_ad > *const _gradp_qp
 Grad(p) at the quadpoints. More...
 
MaterialProperty< std::vector< std::vector< Real > > > *const _dgradp_qp_dgradv
 d(grad porepressure)/d(grad PorousFlow variable) at the quadpoints More...
 
MaterialProperty< std::vector< std::vector< RealGradient > > > *const _dgradp_qp_dv
 d(grad porepressure)/d(PorousFlow variable) at the quadpoints More...
 
GenericMaterialProperty< std::vector< Real >, is_ad > & _saturation
 Computed nodal or qp saturation of the phases. More...
 
MaterialProperty< std::vector< std::vector< Real > > > *const _dsaturation_dvar
 d(saturation)/d(PorousFlow variable) More...
 
GenericMaterialProperty< std::vector< RealGradient >, is_ad > *const _grads_qp
 Grad(s) at the quadpoints. More...
 
MaterialProperty< std::vector< std::vector< Real > > > *const _dgrads_qp_dgradv
 d(grad saturation)/d(grad PorousFlow variable) at the quadpoints More...
 
MaterialProperty< std::vector< std::vector< RealGradient > > > *const _dgrads_qp_dv
 d(grad saturation)/d(PorousFlow variable) at the quadpoints More...
 

Detailed Description

template<bool is_ad>
class PorousFlowFluidStateTempl< is_ad >

Fluid state class using a persistent set of primary variables for the mutliphase, multicomponent case.

Primary variables are: gas pressure, total mass fraction of a component summed over all phases (and optionally temperature in a non-isothermal case).

The total mass fraction of component i summed over all phases, Z_i, is defined as (for two phases)

Z_i = (S_g rho_g Y_i + S_l rho_l X_i) / (S_g rho_g + S_l rho_l)

where S is saturation, rho is density, and the subscripts correspond to gas and liquid phases, respectively, and Y_i and X_i are the mass fractions of the ith component in the gas and liquid phase, respectively.

Depending on the phase conditions, the primary variable Z_i can represent either a mass fraction (when only a single phase is present), or a saturation when two phases are present, and hence it is a persistent variable.

The PorousFlow kernels expect saturation and mass fractions (as well as pressure and temperature), so these must be calculated from Z_i once the state of the system is determined.

A compositional flash calculation using the Rachford-Rice equation is solved to determine vapor fraction (gas saturation), and subsequently the composition of each phase.

Definition at line 47 of file PorousFlowFluidState.h.

Constructor & Destructor Documentation

◆ PorousFlowFluidStateTempl()

template<bool is_ad>
PorousFlowFluidStateTempl< is_ad >::PorousFlowFluidStateTempl ( const InputParameters parameters)

Definition at line 33 of file PorousFlowFluidState.C.

35  _gas_porepressure(_nodal_material
36  ? this->template coupledGenericDofValue<is_ad>("gas_porepressure")
37  : this->template coupledGenericValue<is_ad>("gas_porepressure")),
38  _gas_gradp_qp(this->template coupledGenericGradient<is_ad>("gas_porepressure")),
39  _gas_porepressure_varnum(coupled("gas_porepressure")),
40  _pvar(_dictator.isPorousFlowVariable(_gas_porepressure_varnum)
41  ? _dictator.porousFlowVariableNum(_gas_porepressure_varnum)
42  : 0),
43  _num_Z_vars(coupledComponents("z")),
44  _is_Xnacl_nodal(isCoupled("xnacl") ? getFieldVar("xnacl", 0)->isNodal() : false),
45  _Xnacl(_nodal_material && _is_Xnacl_nodal
46  ? this->template coupledGenericDofValue<is_ad>("xnacl")
47  : this->template coupledGenericValue<is_ad>("xnacl")),
48  _grad_Xnacl_qp(this->template coupledGenericGradient<is_ad>("xnacl")),
49  _Xnacl_varnum(coupled("xnacl")),
50  _Xvar(_dictator.isPorousFlowVariable(_Xnacl_varnum)
51  ? _dictator.porousFlowVariableNum(_Xnacl_varnum)
52  : 0),
53  _fs(this->template getUserObject<PorousFlowFluidStateMultiComponentBase>("fluid_state")),
54  _aqueous_phase_number(_fs.aqueousPhaseIndex()),
55  _gas_phase_number(_fs.gasPhaseIndex()),
56  _aqueous_fluid_component(_fs.aqueousComponentIndex()),
57  _gas_fluid_component(_fs.gasComponentIndex()),
58  _salt_component(_fs.saltComponentIndex()),
60  this->template getGenericMaterialProperty<Real, is_ad>("PorousFlow_temperature" + _sfx)),
61  _gradT_qp(_nodal_material ? nullptr
62  : &this->template getGenericMaterialProperty<RealGradient, is_ad>(
63  "PorousFlow_grad_temperature" + _sfx)),
64  _dtemperature_dvar(is_ad ? nullptr
65  : &this->template getMaterialProperty<std::vector<Real>>(
66  "dPorousFlow_temperature" + _sfx + "_dvar")),
67  _temperature_varnum(coupled("temperature")),
68  _Tvar(_dictator.isPorousFlowVariable(_temperature_varnum)
69  ? _dictator.porousFlowVariableNum(_temperature_varnum)
70  : 0),
71  _pidx(_fs.getPressureIndex()),
72  _Tidx(_fs.getTemperatureIndex()),
73  _Zidx(_fs.getZIndex()),
74  _Xidx(_fs.getXIndex())
75 {
76  // Check that the number of phases in the fluidstate class is also provided in the Dictator
77  if (_fs.numPhases() != _num_phases)
78  mooseError(name(),
79  ": only ",
80  _fs.numPhases(),
81  " phases are allowed. Please check the number of phases entered in the dictator is "
82  "correct");
83 
84  // Store all total mass fractions and associated variable numbers
85  _Z.resize(_num_Z_vars);
86  _gradZ_qp.resize(_num_Z_vars);
87  _Z_varnum.resize(_num_Z_vars);
88  _Zvar.resize(_num_Z_vars);
89 
90  for (unsigned int i = 0; i < _num_Z_vars; ++i)
91  {
92  _Z[i] = (_nodal_material ? &this->template coupledGenericDofValue<is_ad>("z", i)
93  : &this->template coupledGenericValue<is_ad>("z", i));
94  _gradZ_qp[i] = &this->template coupledGenericGradient<is_ad>("z", i);
95  _Z_varnum[i] = coupled("z", i);
96  _Zvar[i] = (_dictator.isPorousFlowVariable(_Z_varnum[i])
97  ? _dictator.porousFlowVariableNum(_Z_varnum[i])
98  : 0);
99  }
100 }
std::vector< unsigned int > _Zvar
PorousFlow variable number of Z.
const GenericVariableValue< is_ad > & _gas_porepressure
Porepressure.
const unsigned int _pidx
Index of derivative wrt pressure.
void mooseError(Args &&... args)
const unsigned int _Tvar
PorousFlow variable number of the temperature.
const MaterialProperty< std::vector< Real > > *const _dtemperature_dvar
Derivative of temperature wrt PorousFlow variables.
const unsigned int _gas_fluid_component
Fluid component number of the gas phase.
const unsigned int _salt_component
Salt component index.
const unsigned int _Xidx
Index of derivative wrt salt mass fraction X.
Compositional flash routines for miscible multiphase flow classes with multiple fluid components...
const unsigned int _aqueous_fluid_component
Fluid component number of the aqueous component.
const unsigned int _Zidx
Index of derivative wrt total mass fraction Z.
const unsigned int _Tidx
Index of derivative wrt temperature.
std::vector< const GenericVariableGradient< is_ad > * > _gradZ_qp
Gradient(s) of total mass fraction(s) of the gas component(s) (only defined at the qps) ...
const std::string _sfx
Suffix to append to material property names (either qp or nodal as required)
const PorousFlowFluidStateMultiComponentBase & _fs
FluidState UserObject.
const std::string name
Definition: Setup.h:20
std::vector< unsigned int > _Z_varnum
Moose variable number of Z.
const unsigned int _gas_phase_number
Phase number of the gas phase.
const unsigned int _gas_porepressure_varnum
Moose variable number of the gas porepressure.
const unsigned int _pvar
PorousFlow variable number of the gas porepressure.
std::vector< const GenericVariableValue< is_ad > * > _Z
Total mass fraction(s) of the gas component(s) summed over all phases.
const GenericVariableValue< is_ad > & _Xnacl
Salt mass fraction (kg/kg)
Fluid state base class using a persistent set of primary variables for multiphase, single and multicomponent cases.
const unsigned int _num_phases
Number of phases.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const GenericVariableGradient< is_ad > & _gas_gradp_qp
Gradient of porepressure (only defined at the qps)
unsigned int numPhases() const
The maximum number of phases in this model.
const unsigned int _num_Z_vars
Number of coupled total mass fractions. Should be _num_phases - 1.
const unsigned int _aqueous_phase_number
Phase number of the aqueous phase.
const unsigned int _temperature_varnum
Moose variable number of the temperature.
const unsigned int _Xnacl_varnum
Salt mass fraction variable number.
const unsigned int _Xvar
Salt mass fraction PorousFlow variable number.
const GenericVariableGradient< is_ad > & _grad_Xnacl_qp
Gradient of salt mass fraction (only defined at the qps)
const GenericMaterialProperty< Real, is_ad > & _temperature
Temperature.
const GenericMaterialProperty< RealGradient, is_ad > *const _gradT_qp
Gradient of temperature (only defined at the qps)
const bool _is_Xnacl_nodal
Flag for nodal NaCl mass fraction.

Member Function Documentation

◆ computeQpProperties()

template<bool is_ad>
void PorousFlowFluidStateTempl< is_ad >::computeQpProperties ( )
overrideprotectedvirtual

Reimplemented from PorousFlowFluidStateBaseMaterialTempl< is_ad >.

Definition at line 122 of file PorousFlowFluidState.C.

123 {
125 
126  // If the material isn't AD, we need to compute the derivatives
127  if (!is_ad)
128  {
129  // Derivative of properties wrt variables (calculated in fluid state class)
130  for (unsigned int ph = 0; ph < _num_phases; ++ph)
131  {
132  // If porepressure is a PorousFlow variable (it usually is), add derivatives wrt
133  // porepressure
134  if (_dictator.isPorousFlowVariable(_gas_porepressure_varnum))
135  {
136  (*_dporepressure_dvar)[_qp][ph][_pvar] = _fsp[ph].pressure.derivatives()[_pidx];
137  (*_dsaturation_dvar)[_qp][ph][_pvar] = _fsp[ph].saturation.derivatives()[_pidx];
138  (*_dfluid_density_dvar)[_qp][ph][_pvar] = _fsp[ph].density.derivatives()[_pidx];
139  (*_dfluid_viscosity_dvar)[_qp][ph][_pvar] = _fsp[ph].viscosity.derivatives()[_pidx];
140  (*_dfluid_enthalpy_dvar)[_qp][ph][_pvar] = _fsp[ph].enthalpy.derivatives()[_pidx];
141  (*_dfluid_internal_energy_dvar)[_qp][ph][_pvar] =
142  _fsp[ph].internal_energy.derivatives()[_pidx];
143 
144  for (unsigned int comp = 0; comp < _num_components; ++comp)
145  (*_dmass_frac_dvar)[_qp][ph][comp][_pvar] =
146  _fsp[ph].mass_fraction[comp].derivatives()[_pidx];
147  }
148 
149  // If Z is a PorousFlow variable (it usually is), add derivatives wrt Z
150  if (_dictator.isPorousFlowVariable(_Z_varnum[0]))
151  {
152  (*_dporepressure_dvar)[_qp][ph][_Zvar[0]] = _fsp[ph].pressure.derivatives()[_Zidx];
153  (*_dsaturation_dvar)[_qp][ph][_Zvar[0]] = _fsp[ph].saturation.derivatives()[_Zidx];
154  (*_dfluid_density_dvar)[_qp][ph][_Zvar[0]] = _fsp[ph].density.derivatives()[_Zidx];
155  (*_dfluid_viscosity_dvar)[_qp][ph][_Zvar[0]] = _fsp[ph].viscosity.derivatives()[_Zidx];
156  (*_dfluid_enthalpy_dvar)[_qp][ph][_Zvar[0]] = _fsp[ph].enthalpy.derivatives()[_Zidx];
157  (*_dfluid_internal_energy_dvar)[_qp][ph][_Zvar[0]] =
158  _fsp[ph].internal_energy.derivatives()[_Zidx];
159 
160  for (unsigned int comp = 0; comp < _num_components; ++comp)
161  (*_dmass_frac_dvar)[_qp][ph][comp][_Zvar[0]] =
162  _fsp[ph].mass_fraction[comp].derivatives()[_Zidx];
163  }
164 
165  // If temperature is a PorousFlow variable (nonisothermal case), add derivatives wrt
166  // temperature
167  if (_dictator.isPorousFlowVariable(_temperature_varnum))
168  {
169  (*_dporepressure_dvar)[_qp][ph][_Tvar] = _fsp[ph].pressure.derivatives()[_Tidx];
170  (*_dsaturation_dvar)[_qp][ph][_Tvar] = _fsp[ph].saturation.derivatives()[_Tidx];
171  (*_dfluid_density_dvar)[_qp][ph][_Tvar] = _fsp[ph].density.derivatives()[_Tidx];
172  (*_dfluid_viscosity_dvar)[_qp][ph][_Tvar] = _fsp[ph].viscosity.derivatives()[_Tidx];
173  (*_dfluid_enthalpy_dvar)[_qp][ph][_Tvar] = _fsp[ph].enthalpy.derivatives()[_Tidx];
174  (*_dfluid_internal_energy_dvar)[_qp][ph][_Tvar] =
175  _fsp[ph].internal_energy.derivatives()[_Tidx];
176 
177  for (unsigned int comp = 0; comp < _num_components; ++comp)
178  (*_dmass_frac_dvar)[_qp][ph][comp][_Tvar] =
179  _fsp[ph].mass_fraction[comp].derivatives()[_Tidx];
180  }
181 
182  // If Xnacl is a PorousFlow variable, add derivatives wrt Xnacl
183  if (_dictator.isPorousFlowVariable(_Xnacl_varnum))
184  {
185  (*_dporepressure_dvar)[_qp][ph][_Xvar] = _fsp[ph].pressure.derivatives()[_Xidx];
186  (*_dsaturation_dvar)[_qp][ph][_Xvar] = _fsp[ph].saturation.derivatives()[_Xidx];
187  (*_dfluid_density_dvar)[_qp][ph][_Xvar] += _fsp[ph].density.derivatives()[_Xidx];
188  (*_dfluid_viscosity_dvar)[_qp][ph][_Xvar] += _fsp[ph].viscosity.derivatives()[_Xidx];
189  (*_dfluid_enthalpy_dvar)[_qp][ph][_Xvar] = _fsp[ph].enthalpy.derivatives()[_Xidx];
190  (*_dfluid_internal_energy_dvar)[_qp][ph][_Xvar] =
191  _fsp[ph].internal_energy.derivatives()[_Xidx];
192 
193  for (unsigned int comp = 0; comp < _num_components; ++comp)
194  (*_dmass_frac_dvar)[_qp][ph][comp][_Xvar] =
195  _fsp[ph].mass_fraction[comp].derivatives()[_Xidx];
196  }
197  }
198  }
199 
200  // If the material properties are being evaluated at the qps, calculate the gradients as well
201  // Note: only nodal properties are evaluated in initQpStatefulProperties(), so no need to check
202  // _is_initqp flag for qp properties
203  if (!_nodal_material)
204  if constexpr (!is_ad)
205  {
206  // Derivatives of capillary pressure
207  const Real dpc = _pc.dCapillaryPressure(_fsp[_aqueous_phase_number].saturation.value(), _qp);
208  const Real d2pc =
209  _pc.d2CapillaryPressure(_fsp[_aqueous_phase_number].saturation.value(), _qp);
210 
211  // Gradients of saturation and porepressure in all phases
212  (*_grads_qp)[_qp][_gas_phase_number] =
213  (*_dsaturation_dvar)[_qp][_gas_phase_number][_pvar] * _gas_gradp_qp[_qp] +
214  (*_dsaturation_dvar)[_qp][_gas_phase_number][_Zvar[0]] * (*_gradZ_qp[0])[_qp] +
215  (*_dsaturation_dvar)[_qp][_gas_phase_number][_Tvar] * (*_gradT_qp)[_qp];
216  (*_grads_qp)[_qp][_aqueous_phase_number] = -(*_grads_qp)[_qp][_gas_phase_number];
217 
218  (*_gradp_qp)[_qp][_gas_phase_number] = _gas_gradp_qp[_qp];
219  (*_gradp_qp)[_qp][_aqueous_phase_number] =
220  _gas_gradp_qp[_qp] - dpc * (*_grads_qp)[_qp][_aqueous_phase_number];
221 
222  // Gradients of mass fractions for each component in each phase
223  (*_grad_mass_frac_qp)[_qp][_aqueous_phase_number][_aqueous_fluid_component] =
224  _fsp[_aqueous_phase_number].mass_fraction[_aqueous_fluid_component].derivatives()[_pidx] *
225  _gas_gradp_qp[_qp] +
226  _fsp[_aqueous_phase_number].mass_fraction[_aqueous_fluid_component].derivatives()[_Zidx] *
227  (*_gradZ_qp[0])[_qp] +
228  _fsp[_aqueous_phase_number].mass_fraction[_aqueous_fluid_component].derivatives()[_Tidx] *
229  (*_gradT_qp)[_qp];
230  (*_grad_mass_frac_qp)[_qp][_aqueous_phase_number][_gas_fluid_component] =
231  -(*_grad_mass_frac_qp)[_qp][_aqueous_phase_number][_aqueous_fluid_component];
232 
233  (*_grad_mass_frac_qp)[_qp][_gas_phase_number][_aqueous_fluid_component] =
234  _fsp[_gas_phase_number].mass_fraction[_aqueous_fluid_component].derivatives()[_pidx] *
235  _gas_gradp_qp[_qp] +
236  _fsp[_gas_phase_number].mass_fraction[_aqueous_fluid_component].derivatives()[_Zidx] *
237  (*_gradZ_qp[0])[_qp] +
238  _fsp[_gas_phase_number].mass_fraction[_aqueous_fluid_component].derivatives()[_Tidx] *
239  (*_gradT_qp)[_qp];
240  (*_grad_mass_frac_qp)[_qp][_gas_phase_number][_gas_fluid_component] =
241  -(*_grad_mass_frac_qp)[_qp][_gas_phase_number][_aqueous_fluid_component];
242 
243  // Derivatives of gradients wrt variables
244  if (_dictator.isPorousFlowVariable(_gas_porepressure_varnum))
245  {
246  for (unsigned int ph = 0; ph < _num_phases; ++ph)
247  (*_dgradp_qp_dgradv)[_qp][ph][_pvar] = 1.0;
248 
249  (*_dgradp_qp_dgradv)[_qp][_aqueous_phase_number][_pvar] +=
250  -dpc * (*_dsaturation_dvar)[_qp][_aqueous_phase_number][_pvar];
251 
252  (*_dgradp_qp_dv)[_qp][_aqueous_phase_number][_pvar] =
253  -d2pc * (*_grads_qp)[_qp][_aqueous_phase_number] *
254  (*_dsaturation_dvar)[_qp][_aqueous_phase_number][_pvar];
255  }
256 
257  if (_dictator.isPorousFlowVariable(_Z_varnum[0]))
258  {
259  (*_dgradp_qp_dgradv)[_qp][_aqueous_phase_number][_Zvar[0]] =
260  -dpc * (*_dsaturation_dvar)[_qp][_aqueous_phase_number][_Zvar[0]];
261 
262  (*_dgradp_qp_dv)[_qp][_aqueous_phase_number][_Zvar[0]] =
263  -d2pc * (*_grads_qp)[_qp][_aqueous_phase_number] *
264  (*_dsaturation_dvar)[_qp][_aqueous_phase_number][_Zvar[0]];
265  }
266 
267  if (_dictator.isPorousFlowVariable(_temperature_varnum))
268  {
269  (*_dgradp_qp_dgradv)[_qp][_aqueous_phase_number][_Tvar] =
270  -dpc * (*_dsaturation_dvar)[_qp][_aqueous_phase_number][_Tvar];
271 
272  (*_dgradp_qp_dv)[_qp][_aqueous_phase_number][_Tvar] =
273  -d2pc * (*_grads_qp)[_qp][_aqueous_phase_number] *
274  (*_dsaturation_dvar)[_qp][_aqueous_phase_number][_Tvar];
275  }
276 
277  // If Xnacl is a PorousFlow variable, add gradients and derivatives wrt Xnacl
278  if (_dictator.isPorousFlowVariable(_Xnacl_varnum))
279  {
280  (*_dgradp_qp_dgradv)[_qp][_aqueous_phase_number][_Xvar] =
281  -dpc * (*_dsaturation_dvar)[_qp][_aqueous_phase_number][_Xvar];
282 
283  (*_grads_qp)[_qp][_aqueous_phase_number] +=
284  (*_dsaturation_dvar)[_qp][_aqueous_phase_number][_Xvar] * _grad_Xnacl_qp[_qp];
285 
286  (*_grads_qp)[_qp][_gas_phase_number] -=
287  (*_dsaturation_dvar)[_qp][_aqueous_phase_number][_Xvar] * _grad_Xnacl_qp[_qp];
288 
289  (*_gradp_qp)[_qp][_aqueous_phase_number] -=
290  dpc * (*_dsaturation_dvar)[_qp][_aqueous_phase_number][_Xvar] * _grad_Xnacl_qp[_qp];
291 
292  (*_dgradp_qp_dv)[_qp][_aqueous_phase_number][_Xvar] =
293  -d2pc * (*_grads_qp)[_qp][_aqueous_phase_number] *
294  (*_dsaturation_dvar)[_qp][_aqueous_phase_number][_Xvar];
295 
296  (*_grad_mass_frac_qp)[_qp][_aqueous_phase_number][_salt_component] = _grad_Xnacl_qp[_qp];
297  (*_grad_mass_frac_qp)[_qp][_aqueous_phase_number][_aqueous_fluid_component] +=
299  .mass_fraction[_aqueous_fluid_component]
300  .derivatives()[_Xidx] *
301  _grad_Xnacl_qp[_qp];
302  (*_grad_mass_frac_qp)[_qp][_aqueous_phase_number][_gas_fluid_component] -=
304  .mass_fraction[_aqueous_fluid_component]
305  .derivatives()[_Xidx] *
306  _grad_Xnacl_qp[_qp];
307  (*_grad_mass_frac_qp)[_qp][_gas_phase_number][_aqueous_fluid_component] +=
308  _fsp[_gas_phase_number].mass_fraction[_aqueous_fluid_component].derivatives()[_Xidx] *
309  _grad_Xnacl_qp[_qp];
310  (*_grad_mass_frac_qp)[_qp][_gas_phase_number][_gas_fluid_component] -=
311  _fsp[_gas_phase_number].mass_fraction[_aqueous_fluid_component].derivatives()[_Xidx] *
312  _grad_Xnacl_qp[_qp];
313  }
314  }
315 }
std::vector< unsigned int > _Zvar
PorousFlow variable number of Z.
const PorousFlowCapillaryPressure & _pc
Capillary pressure UserObject.
const unsigned int _pidx
Index of derivative wrt pressure.
const unsigned int _Tvar
PorousFlow variable number of the temperature.
MaterialProperty< std::vector< std::vector< Real > > > *const _dgradp_qp_dgradv
d(grad porepressure)/d(grad PorousFlow variable) at the quadpoints
const unsigned int _gas_fluid_component
Fluid component number of the gas phase.
const unsigned int _salt_component
Salt component index.
const unsigned int _Xidx
Index of derivative wrt salt mass fraction X.
const unsigned int _aqueous_fluid_component
Fluid component number of the aqueous component.
const unsigned int _Zidx
Index of derivative wrt total mass fraction Z.
const unsigned int _Tidx
Index of derivative wrt temperature.
std::vector< const GenericVariableGradient< is_ad > * > _gradZ_qp
Gradient(s) of total mass fraction(s) of the gas component(s) (only defined at the qps) ...
MaterialProperty< std::vector< std::vector< std::vector< Real > > > > *const _dmass_frac_dvar
Derivative of the mass fraction matrix with respect to the Porous Flow variables. ...
std::vector< unsigned int > _Z_varnum
Moose variable number of Z.
const unsigned int _gas_phase_number
Phase number of the gas phase.
const unsigned int _gas_porepressure_varnum
Moose variable number of the gas porepressure.
const unsigned int _num_components
Number of components.
const unsigned int _pvar
PorousFlow variable number of the gas porepressure.
const unsigned int _num_phases
Number of phases.
virtual Real d2CapillaryPressure(Real saturation, unsigned qp=0) const
Second derivative of capillary pressure wrt true saturation.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real dCapillaryPressure(Real saturation, unsigned qp=0) const
Derivative of capillary pressure wrt true saturation.
const GenericVariableGradient< is_ad > & _gas_gradp_qp
Gradient of porepressure (only defined at the qps)
std::vector< FluidStateProperties > _fsp
FluidStateProperties data for each fluid component in each fluid phase.
const unsigned int _aqueous_phase_number
Phase number of the aqueous phase.
const unsigned int _temperature_varnum
Moose variable number of the temperature.
const unsigned int _Xnacl_varnum
Salt mass fraction variable number.
const unsigned int _Xvar
Salt mass fraction PorousFlow variable number.
const GenericVariableGradient< is_ad > & _grad_Xnacl_qp
Gradient of salt mass fraction (only defined at the qps)
MaterialProperty< std::vector< std::vector< Real > > > *const _dsaturation_dvar
d(saturation)/d(PorousFlow variable)

◆ genericValue() [1/3]

template<bool is_ad>
GenericReal<is_ad> PorousFlowFluidStateBaseMaterialTempl< is_ad >::genericValue ( const ADReal value)
protectedinherited

Templated method that takes an ADReal value, and returns the raw value when is_ad = false, and the AD value when is_ad = true.

Parameters
valueAD value
Returns
raw or AD value depending on is_ad

◆ genericValue() [2/3]

template<>
GenericReal< false > PorousFlowFluidStateBaseMaterialTempl< false >::genericValue ( const ADReal value)
protectedinherited

Definition at line 76 of file PorousFlowFluidStateBaseMaterial.C.

77 {
78  return MetaPhysicL::raw_value(value);
79 }
auto raw_value(const Eigen::Map< T > &in)

◆ genericValue() [3/3]

template<>
GenericReal< true > PorousFlowFluidStateBaseMaterialTempl< true >::genericValue ( const ADReal value)
protectedinherited

Definition at line 83 of file PorousFlowFluidStateBaseMaterial.C.

84 {
85  return value;
86 }
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)

◆ initQpStatefulProperties()

template<bool is_ad>
void PorousFlowFluidStateTempl< is_ad >::initQpStatefulProperties ( )
overrideprotectedvirtual

◆ setMaterialVectorSize()

template<bool is_ad>
void PorousFlowFluidStateBaseMaterialTempl< is_ad >::setMaterialVectorSize ( ) const
protectedvirtualinherited

Size material property vectors and initialise with zeros.

Reimplemented in PorousFlowFluidStateSingleComponentTempl< is_ad >.

Definition at line 149 of file PorousFlowFluidStateBaseMaterial.C.

Referenced by PorousFlowFluidStateSingleComponentTempl< is_ad >::setMaterialVectorSize().

150 {
151  _fluid_density[_qp].resize(_num_phases, 0.0);
152  _fluid_viscosity[_qp].resize(_num_phases, 0.0);
153  _fluid_enthalpy[_qp].resize(_num_phases, 0.0);
154  _fluid_internal_energy[_qp].resize(_num_phases, 0.0);
155  _mass_frac[_qp].resize(_num_phases);
156 
157  for (unsigned int ph = 0; ph < _num_phases; ++ph)
158  _mass_frac[_qp][ph].resize(_num_components);
159 
160  // Derivatives and gradients are not required in initQpStatefulProperties
161  if (!_is_initqp)
162  {
163  // The gradient of the mass fractions are needed for qp materials and AD materials
164  if (!_nodal_material || is_ad)
165  {
166  (*_grad_mass_frac_qp)[_qp].resize(_num_phases);
167 
168  for (unsigned int ph = 0; ph < _num_phases; ++ph)
169  (*_grad_mass_frac_qp)[_qp][ph].assign(_num_components, RealGradient());
170  }
171 
172  // No derivatives are required for AD materials
173  if (!is_ad)
174  {
175  (*_dfluid_density_dvar)[_qp].resize(_num_phases);
176  (*_dfluid_viscosity_dvar)[_qp].resize(_num_phases);
177  (*_dfluid_enthalpy_dvar)[_qp].resize(_num_phases);
178  (*_dfluid_internal_energy_dvar)[_qp].resize(_num_phases);
179  (*_dmass_frac_dvar)[_qp].resize(_num_phases);
180 
181  for (unsigned int ph = 0; ph < _num_phases; ++ph)
182  {
183  (*_dfluid_density_dvar)[_qp][ph].assign(_num_pf_vars, 0.0);
184  (*_dfluid_viscosity_dvar)[_qp][ph].assign(_num_pf_vars, 0.0);
185  (*_dfluid_enthalpy_dvar)[_qp][ph].assign(_num_pf_vars, 0.0);
186  (*_dfluid_internal_energy_dvar)[_qp][ph].assign(_num_pf_vars, 0.0);
187  (*_dmass_frac_dvar)[_qp][ph].resize(_num_components);
188 
189  for (unsigned int comp = 0; comp < _num_components; ++comp)
190  (*_dmass_frac_dvar)[_qp][ph][comp].assign(_num_pf_vars, 0.0);
191  }
192  }
193  }
194 }
GenericMaterialProperty< std::vector< Real >, is_ad > & _fluid_viscosity
Viscosity of each phase.
GenericMaterialProperty< std::vector< Real >, is_ad > & _fluid_internal_energy
Internal energy of each phase.
GenericMaterialProperty< std::vector< std::vector< Real > >, is_ad > & _mass_frac
Mass fraction matrix (indexing is fluid component in fluid phase)
GenericMaterialProperty< std::vector< std::vector< RealGradient > >, is_ad > *const _grad_mass_frac_qp
Gradient of the mass fraction matrix (only defined at the qps)
MaterialProperty< std::vector< std::vector< std::vector< Real > > > > *const _dmass_frac_dvar
Derivative of the mass fraction matrix with respect to the Porous Flow variables. ...
const unsigned int _num_components
Number of components.
const unsigned int _num_pf_vars
Number of PorousFlow variables.
GenericMaterialProperty< std::vector< Real >, is_ad > & _fluid_density
Fluid density of each phase.
bool _is_initqp
Flag to indicate whether stateful properties should be computed.
const unsigned int _num_phases
Number of phases.
GenericMaterialProperty< std::vector< Real >, is_ad > & _fluid_enthalpy
Enthalpy of each phase.

◆ thermophysicalProperties()

template<bool is_ad>
void PorousFlowFluidStateTempl< is_ad >::thermophysicalProperties ( )
overrideprotectedvirtual

Calculates all required thermophysical properties and derivatives for each phase and fluid component.

Implements PorousFlowFluidStateBaseMaterialTempl< is_ad >.

Definition at line 104 of file PorousFlowFluidState.C.

105 {
106  // The FluidProperty objects use temperature in K
107  const GenericReal<is_ad> Tk = _temperature[_qp] + _T_c2k;
108 
110  _fs.thermophysicalProperties(_gas_porepressure[_qp], Tk, _Xnacl[_qp], (*_Z[0])[_qp], _qp, _fsp);
111 }
Moose::GenericType< Real, is_ad > GenericReal
const GenericVariableValue< is_ad > & _gas_porepressure
Porepressure.
const Real _T_c2k
Conversion from degrees Celsius to degrees Kelvin.
void clearFluidStateProperties(std::vector< FluidStateProperties > &fsp) const
Clears the contents of the FluidStateProperties data structure.
virtual void thermophysicalProperties(Real pressure, Real temperature, Real Xnacl, Real Z, unsigned int qp, std::vector< FluidStateProperties > &fsp) const =0
Determines the complete thermophysical state of the system for a given set of primary variables...
const PorousFlowFluidStateMultiComponentBase & _fs
FluidState UserObject.
std::vector< const GenericVariableValue< is_ad > * > _Z
Total mass fraction(s) of the gas component(s) summed over all phases.
const GenericVariableValue< is_ad > & _Xnacl
Salt mass fraction (kg/kg)
std::vector< FluidStateProperties > _fsp
FluidStateProperties data for each fluid component in each fluid phase.
const GenericMaterialProperty< Real, is_ad > & _temperature
Temperature.

◆ validParams()

template<bool is_ad>
InputParameters PorousFlowFluidStateTempl< is_ad >::validParams ( )
static

Definition at line 17 of file PorousFlowFluidState.C.

18 {
20  params.addRequiredCoupledVar("gas_porepressure",
21  "Variable that is the porepressure of the gas phase");
22  params.addRequiredCoupledVar("z", "Total mass fraction of component i summed over all phases");
23  params.addCoupledVar(
24  "temperature", 20, "The fluid temperature (C or K, depending on temperature_unit)");
25  params.addCoupledVar("xnacl", 0, "The salt mass fraction in the brine (kg/kg)");
26  params.addRequiredParam<UserObjectName>("fluid_state", "Name of the FluidState UserObject");
27  params.addClassDescription("Class for fluid state calculations using persistent primary "
28  "variables and a vapor-liquid flash");
29  return params;
30 }
void addRequiredParam(const std::string &name, const std::string &doc_string)
void addCoupledVar(const std::string &name, const std::string &doc_string)
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
void addClassDescription(const std::string &doc_string)

Member Data Documentation

◆ _aqueous_fluid_component

template<bool is_ad>
const unsigned int PorousFlowFluidStateTempl< is_ad >::_aqueous_fluid_component
protected

Fluid component number of the aqueous component.

Definition at line 94 of file PorousFlowFluidState.h.

◆ _aqueous_phase_number

template<bool is_ad>
const unsigned int PorousFlowFluidStateTempl< is_ad >::_aqueous_phase_number
protected

Phase number of the aqueous phase.

Definition at line 90 of file PorousFlowFluidState.h.

◆ _dfluid_density_dvar

template<bool is_ad>
MaterialProperty<std::vector<std::vector<Real> > >* const PorousFlowFluidStateBaseMaterialTempl< is_ad >::_dfluid_density_dvar
protectedinherited

Derivative of the fluid density for each phase wrt PorousFlow variables.

Definition at line 70 of file PorousFlowFluidStateBaseMaterial.h.

◆ _dfluid_enthalpy_dvar

template<bool is_ad>
MaterialProperty<std::vector<std::vector<Real> > >* const PorousFlowFluidStateBaseMaterialTempl< is_ad >::_dfluid_enthalpy_dvar
protectedinherited

Derivative of the fluid enthalpy for each phase wrt PorousFlow variables.

Definition at line 78 of file PorousFlowFluidStateBaseMaterial.h.

◆ _dfluid_internal_energy_dvar

template<bool is_ad>
MaterialProperty<std::vector<std::vector<Real> > >* const PorousFlowFluidStateBaseMaterialTempl< is_ad >::_dfluid_internal_energy_dvar
protectedinherited

Derivative of the fluid internal energy for each phase wrt PorousFlow variables.

Definition at line 82 of file PorousFlowFluidStateBaseMaterial.h.

◆ _dfluid_viscosity_dvar

template<bool is_ad>
MaterialProperty<std::vector<std::vector<Real> > >* const PorousFlowFluidStateBaseMaterialTempl< is_ad >::_dfluid_viscosity_dvar
protectedinherited

Derivative of the fluid viscosity for each phase wrt PorousFlow variables.

Definition at line 74 of file PorousFlowFluidStateBaseMaterial.h.

◆ _dgradp_qp_dgradv

template<bool is_ad>
MaterialProperty<std::vector<std::vector<Real> > >* const PorousFlowVariableBaseTempl< is_ad >::_dgradp_qp_dgradv
protectedinherited

d(grad porepressure)/d(grad PorousFlow variable) at the quadpoints

Definition at line 51 of file PorousFlowVariableBase.h.

◆ _dgradp_qp_dv

template<bool is_ad>
MaterialProperty<std::vector<std::vector<RealGradient> > >* const PorousFlowVariableBaseTempl< is_ad >::_dgradp_qp_dv
protectedinherited

d(grad porepressure)/d(PorousFlow variable) at the quadpoints

Definition at line 54 of file PorousFlowVariableBase.h.

◆ _dgrads_qp_dgradv

template<bool is_ad>
MaterialProperty<std::vector<std::vector<Real> > >* const PorousFlowVariableBaseTempl< is_ad >::_dgrads_qp_dgradv
protectedinherited

d(grad saturation)/d(grad PorousFlow variable) at the quadpoints

Definition at line 66 of file PorousFlowVariableBase.h.

◆ _dgrads_qp_dv

template<bool is_ad>
MaterialProperty<std::vector<std::vector<RealGradient> > >* const PorousFlowVariableBaseTempl< is_ad >::_dgrads_qp_dv
protectedinherited

d(grad saturation)/d(PorousFlow variable) at the quadpoints

Definition at line 69 of file PorousFlowVariableBase.h.

◆ _dmass_frac_dvar

template<bool is_ad>
MaterialProperty<std::vector<std::vector<std::vector<Real> > > >* const PorousFlowFluidStateBaseMaterialTempl< is_ad >::_dmass_frac_dvar
protectedinherited

Derivative of the mass fraction matrix with respect to the Porous Flow variables.

Definition at line 66 of file PorousFlowFluidStateBaseMaterial.h.

◆ _dporepressure_dvar

template<bool is_ad>
MaterialProperty<std::vector<std::vector<Real> > >* const PorousFlowVariableBaseTempl< is_ad >::_dporepressure_dvar
protectedinherited

d(porepressure)/d(PorousFlow variable)

Definition at line 45 of file PorousFlowVariableBase.h.

◆ _dsaturation_dvar

template<bool is_ad>
MaterialProperty<std::vector<std::vector<Real> > >* const PorousFlowVariableBaseTempl< is_ad >::_dsaturation_dvar
protectedinherited

d(saturation)/d(PorousFlow variable)

Definition at line 60 of file PorousFlowVariableBase.h.

◆ _dtemperature_dvar

template<bool is_ad>
const MaterialProperty<std::vector<Real> >* const PorousFlowFluidStateTempl< is_ad >::_dtemperature_dvar
protected

Derivative of temperature wrt PorousFlow variables.

Definition at line 104 of file PorousFlowFluidState.h.

◆ _fluid_density

template<bool is_ad>
GenericMaterialProperty<std::vector<Real>, is_ad>& PorousFlowFluidStateBaseMaterialTempl< is_ad >::_fluid_density
protectedinherited

Fluid density of each phase.

Definition at line 68 of file PorousFlowFluidStateBaseMaterial.h.

◆ _fluid_enthalpy

template<bool is_ad>
GenericMaterialProperty<std::vector<Real>, is_ad>& PorousFlowFluidStateBaseMaterialTempl< is_ad >::_fluid_enthalpy
protectedinherited

Enthalpy of each phase.

Definition at line 76 of file PorousFlowFluidStateBaseMaterial.h.

◆ _fluid_internal_energy

template<bool is_ad>
GenericMaterialProperty<std::vector<Real>, is_ad>& PorousFlowFluidStateBaseMaterialTempl< is_ad >::_fluid_internal_energy
protectedinherited

Internal energy of each phase.

Definition at line 80 of file PorousFlowFluidStateBaseMaterial.h.

◆ _fluid_viscosity

template<bool is_ad>
GenericMaterialProperty<std::vector<Real>, is_ad>& PorousFlowFluidStateBaseMaterialTempl< is_ad >::_fluid_viscosity
protectedinherited

Viscosity of each phase.

Definition at line 72 of file PorousFlowFluidStateBaseMaterial.h.

◆ _fs

template<bool is_ad>
const PorousFlowFluidStateMultiComponentBase& PorousFlowFluidStateTempl< is_ad >::_fs
protected

◆ _fsp

template<bool is_ad>
std::vector<FluidStateProperties> PorousFlowFluidStateBaseMaterialTempl< is_ad >::_fsp
protectedinherited

FluidStateProperties data for each fluid component in each fluid phase.

Definition at line 55 of file PorousFlowFluidStateBaseMaterial.h.

Referenced by PorousFlowFluidStateBaseMaterialTempl< is_ad >::PorousFlowFluidStateBaseMaterialTempl().

◆ _gas_fluid_component

template<bool is_ad>
const unsigned int PorousFlowFluidStateTempl< is_ad >::_gas_fluid_component
protected

Fluid component number of the gas phase.

Definition at line 96 of file PorousFlowFluidState.h.

◆ _gas_gradp_qp

template<bool is_ad>
const GenericVariableGradient<is_ad>& PorousFlowFluidStateTempl< is_ad >::_gas_gradp_qp
protected

Gradient of porepressure (only defined at the qps)

Definition at line 62 of file PorousFlowFluidState.h.

◆ _gas_phase_number

template<bool is_ad>
const unsigned int PorousFlowFluidStateTempl< is_ad >::_gas_phase_number
protected

Phase number of the gas phase.

Definition at line 92 of file PorousFlowFluidState.h.

◆ _gas_porepressure

template<bool is_ad>
const GenericVariableValue<is_ad>& PorousFlowFluidStateTempl< is_ad >::_gas_porepressure
protected

Porepressure.

Definition at line 60 of file PorousFlowFluidState.h.

◆ _gas_porepressure_varnum

template<bool is_ad>
const unsigned int PorousFlowFluidStateTempl< is_ad >::_gas_porepressure_varnum
protected

Moose variable number of the gas porepressure.

Definition at line 64 of file PorousFlowFluidState.h.

◆ _grad_mass_frac_qp

template<bool is_ad>
GenericMaterialProperty<std::vector<std::vector<RealGradient> >, is_ad>* const PorousFlowFluidStateBaseMaterialTempl< is_ad >::_grad_mass_frac_qp
protectedinherited

Gradient of the mass fraction matrix (only defined at the qps)

Definition at line 64 of file PorousFlowFluidStateBaseMaterial.h.

◆ _grad_Xnacl_qp

template<bool is_ad>
const GenericVariableGradient<is_ad>& PorousFlowFluidStateTempl< is_ad >::_grad_Xnacl_qp
protected

Gradient of salt mass fraction (only defined at the qps)

Definition at line 82 of file PorousFlowFluidState.h.

◆ _gradp_qp

template<bool is_ad>
GenericMaterialProperty<std::vector<RealGradient>, is_ad>* const PorousFlowVariableBaseTempl< is_ad >::_gradp_qp
protectedinherited

Grad(p) at the quadpoints.

Definition at line 48 of file PorousFlowVariableBase.h.

◆ _grads_qp

template<bool is_ad>
GenericMaterialProperty<std::vector<RealGradient>, is_ad>* const PorousFlowVariableBaseTempl< is_ad >::_grads_qp
protectedinherited

Grad(s) at the quadpoints.

Definition at line 63 of file PorousFlowVariableBase.h.

◆ _gradT_qp

template<bool is_ad>
const GenericMaterialProperty<RealGradient, is_ad>* const PorousFlowFluidStateTempl< is_ad >::_gradT_qp
protected

Gradient of temperature (only defined at the qps)

Definition at line 102 of file PorousFlowFluidState.h.

◆ _gradZ_qp

template<bool is_ad>
std::vector<const GenericVariableGradient<is_ad> *> PorousFlowFluidStateTempl< is_ad >::_gradZ_qp
protected

Gradient(s) of total mass fraction(s) of the gas component(s) (only defined at the qps)

Definition at line 70 of file PorousFlowFluidState.h.

Referenced by PorousFlowFluidStateTempl< is_ad >::PorousFlowFluidStateTempl().

◆ _is_initqp

template<bool is_ad>
bool PorousFlowFluidStateBaseMaterialTempl< is_ad >::_is_initqp
protectedinherited

Flag to indicate whether stateful properties should be computed.

Definition at line 53 of file PorousFlowFluidStateBaseMaterial.h.

◆ _is_Xnacl_nodal

template<bool is_ad>
const bool PorousFlowFluidStateTempl< is_ad >::_is_Xnacl_nodal
protected

Flag for nodal NaCl mass fraction.

Definition at line 78 of file PorousFlowFluidState.h.

◆ _mass_frac

template<bool is_ad>
GenericMaterialProperty<std::vector<std::vector<Real> >, is_ad>& PorousFlowFluidStateBaseMaterialTempl< is_ad >::_mass_frac
protectedinherited

Mass fraction matrix (indexing is fluid component in fluid phase)

Definition at line 62 of file PorousFlowFluidStateBaseMaterial.h.

◆ _num_components

template<bool is_ad>
const unsigned int PorousFlowVariableBaseTempl< is_ad >::_num_components
protectedinherited

◆ _num_pf_vars

template<bool is_ad>
const unsigned int PorousFlowVariableBaseTempl< is_ad >::_num_pf_vars
protectedinherited

Number of PorousFlow variables.

Definition at line 39 of file PorousFlowVariableBase.h.

◆ _num_phases

template<bool is_ad>
const unsigned int PorousFlowVariableBaseTempl< is_ad >::_num_phases
protectedinherited

◆ _num_Z_vars

template<bool is_ad>
const unsigned int PorousFlowFluidStateTempl< is_ad >::_num_Z_vars
protected

Number of coupled total mass fractions. Should be _num_phases - 1.

Definition at line 76 of file PorousFlowFluidState.h.

Referenced by PorousFlowFluidStateTempl< is_ad >::PorousFlowFluidStateTempl().

◆ _pc

template<bool is_ad>
const PorousFlowCapillaryPressure& PorousFlowFluidStateBaseMaterialTempl< is_ad >::_pc
protectedinherited

Capillary pressure UserObject.

Definition at line 57 of file PorousFlowFluidStateBaseMaterial.h.

◆ _pidx

template<bool is_ad>
const unsigned int PorousFlowFluidStateTempl< is_ad >::_pidx
protected

Index of derivative wrt pressure.

Definition at line 110 of file PorousFlowFluidState.h.

◆ _porepressure

template<bool is_ad>
GenericMaterialProperty<std::vector<Real>, is_ad>& PorousFlowVariableBaseTempl< is_ad >::_porepressure
protectedinherited

◆ _pvar

template<bool is_ad>
const unsigned int PorousFlowFluidStateTempl< is_ad >::_pvar
protected

PorousFlow variable number of the gas porepressure.

Definition at line 66 of file PorousFlowFluidState.h.

◆ _salt_component

template<bool is_ad>
const unsigned int PorousFlowFluidStateTempl< is_ad >::_salt_component
protected

Salt component index.

Definition at line 98 of file PorousFlowFluidState.h.

◆ _saturation

template<bool is_ad>
GenericMaterialProperty<std::vector<Real>, is_ad>& PorousFlowVariableBaseTempl< is_ad >::_saturation
protectedinherited

◆ _sfx

template<bool is_ad>
const std::string PorousFlowFluidStateBaseMaterialTempl< is_ad >::_sfx
protectedinherited

Suffix to append to material property names (either qp or nodal as required)

Definition at line 59 of file PorousFlowFluidStateBaseMaterial.h.

◆ _T_c2k

template<bool is_ad>
const Real PorousFlowFluidStateBaseMaterialTempl< is_ad >::_T_c2k
protectedinherited

Conversion from degrees Celsius to degrees Kelvin.

Definition at line 51 of file PorousFlowFluidStateBaseMaterial.h.

◆ _temperature

template<bool is_ad>
const GenericMaterialProperty<Real, is_ad>& PorousFlowFluidStateTempl< is_ad >::_temperature
protected

Temperature.

Definition at line 100 of file PorousFlowFluidState.h.

◆ _temperature_varnum

template<bool is_ad>
const unsigned int PorousFlowFluidStateTempl< is_ad >::_temperature_varnum
protected

Moose variable number of the temperature.

Definition at line 106 of file PorousFlowFluidState.h.

◆ _Tidx

template<bool is_ad>
const unsigned int PorousFlowFluidStateTempl< is_ad >::_Tidx
protected

Index of derivative wrt temperature.

Definition at line 112 of file PorousFlowFluidState.h.

◆ _Tvar

template<bool is_ad>
const unsigned int PorousFlowFluidStateTempl< is_ad >::_Tvar
protected

PorousFlow variable number of the temperature.

Definition at line 108 of file PorousFlowFluidState.h.

◆ _Xidx

template<bool is_ad>
const unsigned int PorousFlowFluidStateTempl< is_ad >::_Xidx
protected

Index of derivative wrt salt mass fraction X.

Definition at line 116 of file PorousFlowFluidState.h.

◆ _Xnacl

template<bool is_ad>
const GenericVariableValue<is_ad>& PorousFlowFluidStateTempl< is_ad >::_Xnacl
protected

Salt mass fraction (kg/kg)

Definition at line 80 of file PorousFlowFluidState.h.

◆ _Xnacl_varnum

template<bool is_ad>
const unsigned int PorousFlowFluidStateTempl< is_ad >::_Xnacl_varnum
protected

Salt mass fraction variable number.

Definition at line 84 of file PorousFlowFluidState.h.

◆ _Xvar

template<bool is_ad>
const unsigned int PorousFlowFluidStateTempl< is_ad >::_Xvar
protected

Salt mass fraction PorousFlow variable number.

Definition at line 86 of file PorousFlowFluidState.h.

◆ _Z

template<bool is_ad>
std::vector<const GenericVariableValue<is_ad> *> PorousFlowFluidStateTempl< is_ad >::_Z
protected

Total mass fraction(s) of the gas component(s) summed over all phases.

Definition at line 68 of file PorousFlowFluidState.h.

Referenced by PorousFlowFluidStateTempl< is_ad >::PorousFlowFluidStateTempl().

◆ _Z_varnum

template<bool is_ad>
std::vector<unsigned int> PorousFlowFluidStateTempl< is_ad >::_Z_varnum
protected

Moose variable number of Z.

Definition at line 72 of file PorousFlowFluidState.h.

Referenced by PorousFlowFluidStateTempl< is_ad >::PorousFlowFluidStateTempl().

◆ _Zidx

template<bool is_ad>
const unsigned int PorousFlowFluidStateTempl< is_ad >::_Zidx
protected

Index of derivative wrt total mass fraction Z.

Definition at line 114 of file PorousFlowFluidState.h.

◆ _Zvar

template<bool is_ad>
std::vector<unsigned int> PorousFlowFluidStateTempl< is_ad >::_Zvar
protected

PorousFlow variable number of Z.

Definition at line 74 of file PorousFlowFluidState.h.

Referenced by PorousFlowFluidStateTempl< is_ad >::PorousFlowFluidStateTempl().

◆ usingPorousFlowFluidStateBaseMaterialDerivativeMembers

template<bool is_ad>
PorousFlowFluidStateTempl< is_ad >::usingPorousFlowFluidStateBaseMaterialDerivativeMembers
protected

Definition at line 121 of file PorousFlowFluidState.h.

◆ usingPorousFlowFluidStateBaseMaterialMembers

template<bool is_ad>
PorousFlowFluidStateTempl< is_ad >::usingPorousFlowFluidStateBaseMaterialMembers
protected

Definition at line 119 of file PorousFlowFluidState.h.

◆ usingPorousFlowVariableBaseMembers

template<bool is_ad>
PorousFlowFluidStateBaseMaterialTempl< is_ad >::usingPorousFlowVariableBaseMembers
protectedinherited

Definition at line 84 of file PorousFlowFluidStateBaseMaterial.h.


The documentation for this class was generated from the following files: