www.mooseframework.org
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
PorousFlowFluidStateSingleComponent Class Reference

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

#include <PorousFlowFluidStateSingleComponent.h>

Inheritance diagram for PorousFlowFluidStateSingleComponent:
[legend]

Public Member Functions

 PorousFlowFluidStateSingleComponent (const InputParameters &parameters)
 

Protected Member Functions

virtual void initQpStatefulProperties () override
 
virtual void computeQpProperties () override
 
void setMaterialVectorSize () const
 Size material property vectors and initialise with zeros. More...
 
virtual void thermophysicalProperties ()
 Calculates all required thermophysical properties and derivatives for each phase and fluid component. More...
 

Protected Attributes

const VariableValue & _liquid_porepressure
 Porepressure. More...
 
const VariableGradient & _liquid_gradp_qp
 Gradient of porepressure (only defined at the qps) More...
 
const unsigned int _liquid_porepressure_varnum
 Moose variable number of the porepressure. More...
 
const unsigned int _pvar
 PorousFlow variable number of the porepressure. More...
 
const VariableValue & _enthalpy
 Enthalpy. More...
 
const VariableGradient & _gradh_qp
 Gradient of enthalpy (only defined at the qps) More...
 
const unsigned int _enthalpy_varnum
 Moose variable number of the enthalpy. More...
 
const unsigned int _hvar
 PorousFlow variable number of the enthalpy. More...
 
const PorousFlowFluidStateSingleComponentBase_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...
 
MaterialProperty< Real > & _temperature
 Temperature. More...
 
MaterialProperty< RealGradient > * _grad_temperature_qp
 Gradient of temperature (only defined at the qps) More...
 
MaterialProperty< std::vector< Real > > & _dtemperature_dvar
 Derivative of temperature wrt PorousFlow variables. More...
 
MaterialProperty< std::vector< Real > > *const _dgrad_temperature_dgradv
 d(grad temperature)/d(grad PorousFlow variable) at the quadpoints More...
 
MaterialProperty< std::vector< RealGradient > > *const _dgrad_temperature_dv
 d(grad temperature)/d(PorousFlow variable) at the quadpoints More...
 
MaterialProperty< std::vector< std::vector< Real > > > & _mass_frac
 Mass fraction matrix. More...
 
MaterialProperty< std::vector< std::vector< RealGradient > > > * _grad_mass_frac_qp
 Gradient of the mass fraction matrix (only defined at the qps) More...
 
MaterialProperty< std::vector< std::vector< std::vector< Real > > > > & _dmass_frac_dvar
 Derivative of the mass fraction matrix with respect to the Porous Flow variables. More...
 
MaterialProperty< std::vector< Real > > & _fluid_density
 Fluid density of each phase. More...
 
MaterialProperty< std::vector< std::vector< Real > > > & _dfluid_density_dvar
 Derivative of the fluid density for each phase wrt PorousFlow variables. More...
 
MaterialProperty< std::vector< Real > > & _fluid_viscosity
 Viscosity of each phase. More...
 
MaterialProperty< std::vector< std::vector< Real > > > & _dfluid_viscosity_dvar
 Derivative of the fluid viscosity for each phase wrt PorousFlow variables. More...
 
MaterialProperty< std::vector< Real > > & _fluid_enthalpy
 Enthalpy of each phase. More...
 
MaterialProperty< std::vector< std::vector< Real > > > & _dfluid_enthalpy_dvar
 Derivative of the fluid enthalpy for each phase wrt PorousFlow variables. More...
 
MaterialProperty< std::vector< Real > > & _fluid_internal_energy
 Internal energy of each phase. More...
 
MaterialProperty< std::vector< std::vector< Real > > > & _dfluid_internal_energy_dvar
 Derivative of the fluid internal energy for each phase wrt PorousFlow variables. More...
 
const Real _T_c2k
 Conversion from degrees Celsius to degrees Kelvin. More...
 
bool _is_initqp
 Flag to indicate whether to calculate stateful properties. More...
 
std::vector< FluidStateProperties_fsp
 FluidStateProperties data structure. More...
 
FluidStatePhaseEnum _phase_state
 FluidStatePhaseEnum. More...
 
const PorousFlowCapillaryPressure_pc
 Capillary pressure UserObject. More...
 
const unsigned int _pidx
 Index of derivative wrt pressure. More...
 
const unsigned int _hidx
 Index of derivative wrt enthalpy. More...
 
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...
 
MaterialProperty< std::vector< Real > > & _porepressure
 Computed nodal or quadpoint values of porepressure of the phases. More...
 
MaterialProperty< std::vector< std::vector< Real > > > & _dporepressure_dvar
 d(porepressure)/d(PorousFlow variable) More...
 
MaterialProperty< std::vector< RealGradient > > *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...
 
MaterialProperty< std::vector< Real > > & _saturation
 Computed nodal or qp saturation of the phases. More...
 
MaterialProperty< std::vector< std::vector< Real > > > & _dsaturation_dvar
 d(saturation)/d(PorousFlow variable) More...
 
MaterialProperty< std::vector< RealGradient > > *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

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

Primary variables are: liquid pressure and enthalpy.

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

Definition at line 31 of file PorousFlowFluidStateSingleComponent.h.

Constructor & Destructor Documentation

◆ PorousFlowFluidStateSingleComponent()

PorousFlowFluidStateSingleComponent::PorousFlowFluidStateSingleComponent ( const InputParameters &  parameters)

Definition at line 35 of file PorousFlowFluidStateSingleComponent.C.

37  : PorousFlowVariableBase(parameters),
38 
39  _liquid_porepressure(_nodal_material ? coupledDofValues("porepressure")
40  : coupledValue("porepressure")),
41  _liquid_gradp_qp(coupledGradient("porepressure")),
42  _liquid_porepressure_varnum(coupled("porepressure")),
43  _pvar(_dictator.isPorousFlowVariable(_liquid_porepressure_varnum)
44  ? _dictator.porousFlowVariableNum(_liquid_porepressure_varnum)
45  : 0),
46 
47  _enthalpy(_nodal_material ? coupledDofValues("enthalpy") : coupledValue("enthalpy")),
48  _gradh_qp(coupledGradient("enthalpy")),
49  _enthalpy_varnum(coupled("enthalpy")),
50  _hvar(_dictator.isPorousFlowVariable(_enthalpy_varnum)
51  ? _dictator.porousFlowVariableNum(_enthalpy_varnum)
52  : 0),
53 
54  _fs(getUserObject<PorousFlowFluidStateSingleComponentBase>("fluid_state")),
55  _aqueous_phase_number(_fs.aqueousPhaseIndex()),
56  _gas_phase_number(_fs.gasPhaseIndex()),
57 
58  _temperature(_nodal_material ? declareProperty<Real>("PorousFlow_temperature_nodal")
59  : declareProperty<Real>("PorousFlow_temperature_qp")),
60  _grad_temperature_qp(_nodal_material
61  ? nullptr
62  : &declareProperty<RealGradient>("PorousFlow_grad_temperature_qp")),
64  _nodal_material ? declareProperty<std::vector<Real>>("dPorousFlow_temperature_nodal_dvar")
65  : declareProperty<std::vector<Real>>("dPorousFlow_temperature_qp_dvar")),
66  _dgrad_temperature_dgradv(_nodal_material ? nullptr
67  : &declareProperty<std::vector<Real>>(
68  "dPorousFlow_grad_temperature_qp_dgradvar")),
69  _dgrad_temperature_dv(_nodal_material ? nullptr
70  : &declareProperty<std::vector<RealGradient>>(
71  "dPorousFlow_grad_temperature_qp_dvar")),
72  _mass_frac(_nodal_material
73  ? declareProperty<std::vector<std::vector<Real>>>("PorousFlow_mass_frac_nodal")
74  : declareProperty<std::vector<std::vector<Real>>>("PorousFlow_mass_frac_qp")),
75  _grad_mass_frac_qp(_nodal_material ? nullptr
76  : &declareProperty<std::vector<std::vector<RealGradient>>>(
77  "PorousFlow_grad_mass_frac_qp")),
78  _dmass_frac_dvar(_nodal_material ? declareProperty<std::vector<std::vector<std::vector<Real>>>>(
79  "dPorousFlow_mass_frac_nodal_dvar")
80  : declareProperty<std::vector<std::vector<std::vector<Real>>>>(
81  "dPorousFlow_mass_frac_qp_dvar")),
82 
83  _fluid_density(_nodal_material
84  ? declareProperty<std::vector<Real>>("PorousFlow_fluid_phase_density_nodal")
85  : declareProperty<std::vector<Real>>("PorousFlow_fluid_phase_density_qp")),
86  _dfluid_density_dvar(_nodal_material ? declareProperty<std::vector<std::vector<Real>>>(
87  "dPorousFlow_fluid_phase_density_nodal_dvar")
88  : declareProperty<std::vector<std::vector<Real>>>(
89  "dPorousFlow_fluid_phase_density_qp_dvar")),
90  _fluid_viscosity(_nodal_material
91  ? declareProperty<std::vector<Real>>("PorousFlow_viscosity_nodal")
92  : declareProperty<std::vector<Real>>("PorousFlow_viscosity_qp")),
94  _nodal_material
95  ? declareProperty<std::vector<std::vector<Real>>>("dPorousFlow_viscosity_nodal_dvar")
96  : declareProperty<std::vector<std::vector<Real>>>("dPorousFlow_viscosity_qp_dvar")),
97 
99  _nodal_material
100  ? declareProperty<std::vector<Real>>("PorousFlow_fluid_phase_enthalpy_nodal")
101  : declareProperty<std::vector<Real>>("PorousFlow_fluid_phase_enthalpy_qp")),
102  _dfluid_enthalpy_dvar(_nodal_material ? declareProperty<std::vector<std::vector<Real>>>(
103  "dPorousFlow_fluid_phase_enthalpy_nodal_dvar")
104  : declareProperty<std::vector<std::vector<Real>>>(
105  "dPorousFlow_fluid_phase_enthalpy_qp_dvar")),
106 
108  _nodal_material
109  ? declareProperty<std::vector<Real>>("PorousFlow_fluid_phase_internal_energy_nodal")
110  : declareProperty<std::vector<Real>>("PorousFlow_fluid_phase_internal_energy_qp")),
111  _dfluid_internal_energy_dvar(_nodal_material
112  ? declareProperty<std::vector<std::vector<Real>>>(
113  "dPorousFlow_fluid_phase_internal_energy_nodal_dvar")
114  : declareProperty<std::vector<std::vector<Real>>>(
115  "dPorousFlow_fluid_phase_internal_qp_dvar")),
116 
117  _T_c2k(getParam<MooseEnum>("temperature_unit") == 0 ? 0.0 : 273.15),
118  _is_initqp(false),
119  _pc(getUserObject<PorousFlowCapillaryPressure>("capillary_pressure")),
120  _pidx(_fs.getPressureIndex()),
121  _hidx(_fs.getEnthalpyIndex())
122 {
123  // Check that the number of phases in the fluidstate class is also provided in the Dictator
124  if (_fs.numPhases() != _num_phases)
125  mooseError(name(),
126  ": only ",
127  _fs.numPhases(),
128  " phases are allowed. Please check the number of phases entered in the dictator is "
129  "correct");
130 
131  // Set the size of the FluidStateProperties vector
133 }

Member Function Documentation

◆ computeQpProperties()

void PorousFlowFluidStateSingleComponent::computeQpProperties ( )
overrideprotectedvirtual

Reimplemented from PorousFlowVariableBase.

Definition at line 176 of file PorousFlowFluidStateSingleComponent.C.

177 {
178  _is_initqp = false;
179  // Prepare the derivative vectors
181 
182  // Set the size of all other vectors
184 
185  // Calculate all required thermophysical properties
187 
188  // Temperature doesn't depend on fluid phase
189  _temperature[_qp] = _fsp[_aqueous_phase_number].temperature.value() - _T_c2k;
190  _dtemperature_dvar[_qp][_pvar] = _fsp[_aqueous_phase_number].temperature.derivatives()[_pidx];
191  _dtemperature_dvar[_qp][_hvar] = _fsp[_aqueous_phase_number].temperature.derivatives()[_hidx];
192 
193  for (unsigned int ph = 0; ph < _num_phases; ++ph)
194  {
195  _saturation[_qp][ph] = _fsp[ph].saturation.value();
196  _porepressure[_qp][ph] = _fsp[ph].pressure.value();
197  _fluid_density[_qp][ph] = _fsp[ph].density.value();
198  _fluid_viscosity[_qp][ph] = _fsp[ph].viscosity.value();
199  _fluid_enthalpy[_qp][ph] = _fsp[ph].enthalpy.value();
200  _fluid_internal_energy[_qp][ph] = _fsp[ph].internal_energy.value();
201  }
202 
203  // Derivative of pressure, saturation and fluid properties wrt variables
204  for (unsigned int ph = 0; ph < _num_phases; ++ph)
205  {
206  _dporepressure_dvar[_qp][ph][_pvar] = _fsp[ph].pressure.derivatives()[_pidx];
207  _dporepressure_dvar[_qp][ph][_hvar] = _fsp[ph].pressure.derivatives()[_hidx];
208 
209  _dsaturation_dvar[_qp][ph][_pvar] = _fsp[ph].saturation.derivatives()[_pidx];
210  _dsaturation_dvar[_qp][ph][_hvar] = _fsp[ph].saturation.derivatives()[_hidx];
211 
212  _dfluid_density_dvar[_qp][ph][_pvar] = _fsp[ph].density.derivatives()[_pidx];
213  _dfluid_density_dvar[_qp][ph][_hvar] = _fsp[ph].density.derivatives()[_hidx];
214 
215  _dfluid_viscosity_dvar[_qp][ph][_pvar] = _fsp[ph].viscosity.derivatives()[_pidx];
216  _dfluid_viscosity_dvar[_qp][ph][_hvar] = _fsp[ph].viscosity.derivatives()[_hidx];
217 
218  _dfluid_enthalpy_dvar[_qp][ph][_pvar] = _fsp[ph].enthalpy.derivatives()[_pidx];
219  _dfluid_enthalpy_dvar[_qp][ph][_hvar] = _fsp[ph].enthalpy.derivatives()[_hidx];
220 
221  _dfluid_internal_energy_dvar[_qp][ph][_pvar] = _fsp[ph].internal_energy.derivatives()[_pidx];
222  _dfluid_internal_energy_dvar[_qp][ph][_hvar] = _fsp[ph].internal_energy.derivatives()[_hidx];
223  }
224 
225  // If the material properties are being evaluated at the qps, calculate the
226  // gradients as well. Note: only nodal properties are evaluated in
227  // initQpStatefulProperties(), so no need to check _is_initqp flag for qp
228  // properties
229  if (!_nodal_material)
230  {
231  // Need to compute second derivatives of properties wrt variables for some of
232  // the gradient derivatives. Use finite differences for now
233  const Real dp = 1.0e-2;
234  const Real dh = 1.0e-2;
235 
236  std::vector<FluidStateProperties> fsp_dp(_num_phases, FluidStateProperties(_num_components));
238  _liquid_porepressure[_qp] + dp, _enthalpy[_qp], _qp, _phase_state, fsp_dp);
239 
240  std::vector<FluidStateProperties> fsp_dh(_num_phases, FluidStateProperties(_num_components));
242  _liquid_porepressure[_qp], _enthalpy[_qp] + dh, _qp, _phase_state, fsp_dh);
243 
244  // Gradient of temperature (non-zero in all phases)
245  (*_grad_temperature_qp)[_qp] = _dtemperature_dvar[_qp][_pvar] * _liquid_gradp_qp[_qp] +
246  _dtemperature_dvar[_qp][_hvar] * _gradh_qp[_qp];
247  (*_dgrad_temperature_dgradv)[_qp][_pvar] = _dtemperature_dvar[_qp][_pvar];
248  (*_dgrad_temperature_dgradv)[_qp][_hvar] = _dtemperature_dvar[_qp][_hvar];
249 
250  const Real d2T_dp2 = (fsp_dp[_aqueous_phase_number].temperature.derivatives()[_pidx] -
251  _fsp[_aqueous_phase_number].temperature.derivatives()[_pidx]) /
252  dp;
253 
254  const Real d2T_dh2 = (fsp_dh[_aqueous_phase_number].temperature.derivatives()[_hidx] -
255  _fsp[_aqueous_phase_number].temperature.derivatives()[_hidx]) /
256  dh;
257 
258  const Real d2T_dph = (fsp_dp[_aqueous_phase_number].temperature.derivatives()[_hidx] -
259  _fsp[_aqueous_phase_number].temperature.derivatives()[_hidx]) /
260  (2.0 * dp) +
261  (fsp_dh[_aqueous_phase_number].temperature.derivatives()[_pidx] -
262  _fsp[_aqueous_phase_number].temperature.derivatives()[_pidx]) /
263  (2.0 * dh);
264 
265  (*_dgrad_temperature_dv)[_qp][_pvar] =
266  d2T_dp2 * _liquid_gradp_qp[_qp] + d2T_dph * _gradh_qp[_qp];
267  (*_dgrad_temperature_dv)[_qp][_hvar] =
268  d2T_dph * _liquid_gradp_qp[_qp] + d2T_dh2 * _gradh_qp[_qp];
269 
270  // Gradient of saturation and derivatives
271  (*_grads_qp)[_qp][_gas_phase_number] =
274  (*_grads_qp)[_qp][_aqueous_phase_number] = -(*_grads_qp)[_qp][_gas_phase_number];
275 
276  (*_dgrads_qp_dgradv)[_qp][_gas_phase_number][_pvar] =
278  (*_dgrads_qp_dgradv)[_qp][_aqueous_phase_number][_pvar] =
279  -(*_dgrads_qp_dgradv)[_qp][_gas_phase_number][_pvar];
280 
281  (*_dgrads_qp_dgradv)[_qp][_gas_phase_number][_hvar] =
283  (*_dgrads_qp_dgradv)[_qp][_aqueous_phase_number][_hvar] =
284  -(*_dgrads_qp_dgradv)[_qp][_gas_phase_number][_hvar];
285 
286  const Real d2s_dp2 = (fsp_dp[_gas_phase_number].saturation.derivatives()[_pidx] -
287  _fsp[_gas_phase_number].saturation.derivatives()[_pidx]) /
288  dp;
289 
290  const Real d2s_dh2 = (fsp_dh[_gas_phase_number].saturation.derivatives()[_hidx] -
291  _fsp[_gas_phase_number].saturation.derivatives()[_hidx]) /
292  dh;
293 
294  const Real d2s_dph = (fsp_dp[_gas_phase_number].saturation.derivatives()[_hidx] -
295  _fsp[_gas_phase_number].saturation.derivatives()[_hidx]) /
296  (2.0 * dp) +
297  (fsp_dh[_gas_phase_number].saturation.derivatives()[_pidx] -
298  _fsp[_gas_phase_number].saturation.derivatives()[_pidx]) /
299  (2.0 * dh);
300 
301  (*_dgrads_qp_dv)[_qp][_gas_phase_number][_pvar] =
302  d2s_dp2 * _liquid_gradp_qp[_qp] + d2s_dph * _gradh_qp[_qp];
303  (*_dgrads_qp_dv)[_qp][_aqueous_phase_number][_pvar] =
304  -(*_dgrads_qp_dv)[_qp][_gas_phase_number][_pvar];
305 
306  (*_dgrads_qp_dv)[_qp][_gas_phase_number][_hvar] =
307  d2s_dh2 * _gradh_qp[_qp] + d2s_dph * _liquid_gradp_qp[_qp];
308  (*_dgrads_qp_dv)[_qp][_aqueous_phase_number][_hvar] =
309  -(*_dgrads_qp_dv)[_qp][_gas_phase_number][_hvar];
310 
311  // Gradient of porepressure and derivatives
312  // Note: need first and second derivativea of capillary pressure
313  const Real dpc = _pc.dCapillaryPressure(_fsp[_aqueous_phase_number].saturation.value());
314  const Real d2pc = _pc.d2CapillaryPressure(_fsp[_aqueous_phase_number].saturation.value());
315 
316  (*_gradp_qp)[_qp][_aqueous_phase_number] = _liquid_gradp_qp[_qp];
317  (*_gradp_qp)[_qp][_gas_phase_number] =
318  _liquid_gradp_qp[_qp] + dpc * (*_grads_qp)[_qp][_aqueous_phase_number];
319 
320  for (unsigned int ph = 0; ph < _num_phases; ++ph)
321  (*_dgradp_qp_dgradv)[_qp][ph][_pvar] = 1.0;
322 
323  (*_dgradp_qp_dgradv)[_qp][_gas_phase_number][_pvar] +=
324  dpc * (*_dgrads_qp_dgradv)[_qp][_aqueous_phase_number][_pvar];
325  (*_dgradp_qp_dgradv)[_qp][_gas_phase_number][_hvar] =
326  dpc * (*_dgrads_qp_dgradv)[_qp][_aqueous_phase_number][_hvar];
327 
328  (*_dgradp_qp_dv)[_qp][_gas_phase_number][_pvar] =
329  d2pc * (*_grads_qp)[_qp][_aqueous_phase_number] *
331  dpc * (*_dgrads_qp_dv)[_qp][_aqueous_phase_number][_pvar];
332 
333  (*_dgradp_qp_dv)[_qp][_gas_phase_number][_hvar] =
334  d2pc * (*_grads_qp)[_qp][_aqueous_phase_number] *
336  dpc * (*_dgrads_qp_dv)[_qp][_aqueous_phase_number][_hvar];
337  }
338 }

◆ initQpStatefulProperties()

void PorousFlowFluidStateSingleComponent::initQpStatefulProperties ( )
overrideprotectedvirtual

Reimplemented from PorousFlowVariableBase.

Definition at line 143 of file PorousFlowFluidStateSingleComponent.C.

144 {
145  _is_initqp = true;
146  // Set the size of pressure and saturation vectors
148 
149  // Set the size of all other vectors
151 
153 
154  // Set the initial values of the properties at the nodes.
155  // Note: not required for qp materials as no old values at the qps are requested
156  if (_nodal_material)
157  {
159 
160  // Temperature doesn't depend on fluid phase
161  _temperature[_qp] = _fsp[_aqueous_phase_number].temperature.value() - _T_c2k;
162 
163  for (unsigned int ph = 0; ph < _num_phases; ++ph)
164  {
165  _saturation[_qp][ph] = _fsp[ph].saturation.value();
166  _porepressure[_qp][ph] = _fsp[ph].pressure.value();
167  _fluid_density[_qp][ph] = _fsp[ph].density.value();
168  _fluid_viscosity[_qp][ph] = _fsp[ph].viscosity.value();
169  _fluid_enthalpy[_qp][ph] = _fsp[ph].enthalpy.value();
170  _fluid_internal_energy[_qp][ph] = _fsp[ph].internal_energy.value();
171  }
172  }
173 }

◆ setMaterialVectorSize()

void PorousFlowFluidStateSingleComponent::setMaterialVectorSize ( ) const
protected

Size material property vectors and initialise with zeros.

Definition at line 341 of file PorousFlowFluidStateSingleComponent.C.

342 {
343  _fluid_density[_qp].assign(_num_phases, 0.0);
344  _fluid_viscosity[_qp].assign(_num_phases, 0.0);
345  _fluid_enthalpy[_qp].assign(_num_phases, 0.0);
346  _fluid_internal_energy[_qp].assign(_num_phases, 0.0);
347  _mass_frac[_qp].resize(_num_phases);
348 
349  for (unsigned int ph = 0; ph < _num_phases; ++ph)
350  _mass_frac[_qp][ph].assign(_num_components, 1.0);
351 
352  // Derivatives and gradients are not required in initQpStatefulProperties
353  if (!_is_initqp)
354  {
355  _dfluid_density_dvar[_qp].resize(_num_phases);
356  _dfluid_viscosity_dvar[_qp].resize(_num_phases);
357  _dfluid_enthalpy_dvar[_qp].resize(_num_phases);
359  _dmass_frac_dvar[_qp].resize(_num_phases);
360 
361  // Temperature doesn't depend of fluid phase
362  _dtemperature_dvar[_qp].assign(_num_pf_vars, 0.0);
363 
364  if (!_nodal_material)
365  {
366  (*_grad_mass_frac_qp)[_qp].resize(_num_phases);
367  (*_grad_temperature_qp)[_qp] = RealGradient();
368  (*_dgrad_temperature_dgradv)[_qp].assign(_num_pf_vars, 0.0);
369  (*_dgrad_temperature_dv)[_qp].assign(_num_pf_vars, RealGradient());
370  }
371 
372  for (unsigned int ph = 0; ph < _num_phases; ++ph)
373  {
374  _dfluid_density_dvar[_qp][ph].assign(_num_pf_vars, 0.0);
375  _dfluid_viscosity_dvar[_qp][ph].assign(_num_pf_vars, 0.0);
376  _dfluid_enthalpy_dvar[_qp][ph].assign(_num_pf_vars, 0.0);
377  _dfluid_internal_energy_dvar[_qp][ph].assign(_num_pf_vars, 0.0);
378  _dmass_frac_dvar[_qp][ph].resize(_num_components);
379 
380  for (unsigned int comp = 0; comp < _num_components; ++comp)
381  _dmass_frac_dvar[_qp][ph][comp].assign(_num_pf_vars, 0.0);
382 
383  if (!_nodal_material)
384  (*_grad_mass_frac_qp)[_qp][ph].assign(_num_components, RealGradient());
385  }
386  }
387 }

Referenced by computeQpProperties(), and initQpStatefulProperties().

◆ thermophysicalProperties()

void PorousFlowFluidStateSingleComponent::thermophysicalProperties ( )
protectedvirtual

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

Definition at line 136 of file PorousFlowFluidStateSingleComponent.C.

Referenced by computeQpProperties(), and initQpStatefulProperties().

Member Data Documentation

◆ _aqueous_phase_number

const unsigned int PorousFlowFluidStateSingleComponent::_aqueous_phase_number
protected

Phase number of the aqueous phase.

Definition at line 68 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties(), and initQpStatefulProperties().

◆ _dfluid_density_dvar

MaterialProperty<std::vector<std::vector<Real> > >& PorousFlowFluidStateSingleComponent::_dfluid_density_dvar
protected

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

Definition at line 90 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties(), and setMaterialVectorSize().

◆ _dfluid_enthalpy_dvar

MaterialProperty<std::vector<std::vector<Real> > >& PorousFlowFluidStateSingleComponent::_dfluid_enthalpy_dvar
protected

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

Definition at line 98 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties(), and setMaterialVectorSize().

◆ _dfluid_internal_energy_dvar

MaterialProperty<std::vector<std::vector<Real> > >& PorousFlowFluidStateSingleComponent::_dfluid_internal_energy_dvar
protected

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

Definition at line 102 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties(), and setMaterialVectorSize().

◆ _dfluid_viscosity_dvar

MaterialProperty<std::vector<std::vector<Real> > >& PorousFlowFluidStateSingleComponent::_dfluid_viscosity_dvar
protected

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

Definition at line 94 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties(), and setMaterialVectorSize().

◆ _dgrad_temperature_dgradv

MaterialProperty<std::vector<Real> >* const PorousFlowFluidStateSingleComponent::_dgrad_temperature_dgradv
protected

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

Definition at line 78 of file PorousFlowFluidStateSingleComponent.h.

◆ _dgrad_temperature_dv

MaterialProperty<std::vector<RealGradient> >* const PorousFlowFluidStateSingleComponent::_dgrad_temperature_dv
protected

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

Definition at line 80 of file PorousFlowFluidStateSingleComponent.h.

◆ _dgradp_qp_dgradv

MaterialProperty<std::vector<std::vector<Real> > >* const PorousFlowVariableBase::_dgradp_qp_dgradv
protectedinherited

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

Definition at line 53 of file PorousFlowVariableBase.h.

Referenced by computeQpProperties(), and PorousFlowFluidState::computeQpProperties().

◆ _dgradp_qp_dv

MaterialProperty<std::vector<std::vector<RealGradient> > >* const PorousFlowVariableBase::_dgradp_qp_dv
protectedinherited

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

Definition at line 56 of file PorousFlowVariableBase.h.

◆ _dgrads_qp_dgradv

MaterialProperty<std::vector<std::vector<Real> > >* const PorousFlowVariableBase::_dgrads_qp_dgradv
protectedinherited

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

Definition at line 68 of file PorousFlowVariableBase.h.

◆ _dgrads_qp_dv

MaterialProperty<std::vector<std::vector<RealGradient> > >* const PorousFlowVariableBase::_dgrads_qp_dv
protectedinherited

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

Definition at line 71 of file PorousFlowVariableBase.h.

◆ _dmass_frac_dvar

MaterialProperty<std::vector<std::vector<std::vector<Real> > > >& PorousFlowFluidStateSingleComponent::_dmass_frac_dvar
protected

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

Definition at line 86 of file PorousFlowFluidStateSingleComponent.h.

Referenced by setMaterialVectorSize().

◆ _dporepressure_dvar

MaterialProperty<std::vector<std::vector<Real> > >& PorousFlowVariableBase::_dporepressure_dvar
protectedinherited

◆ _dsaturation_dvar

MaterialProperty<std::vector<std::vector<Real> > >& PorousFlowVariableBase::_dsaturation_dvar
protectedinherited

◆ _dtemperature_dvar

MaterialProperty<std::vector<Real> >& PorousFlowFluidStateSingleComponent::_dtemperature_dvar
protected

Derivative of temperature wrt PorousFlow variables.

Definition at line 76 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties(), and setMaterialVectorSize().

◆ _enthalpy

const VariableValue& PorousFlowFluidStateSingleComponent::_enthalpy
protected

Enthalpy.

Definition at line 58 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties(), and thermophysicalProperties().

◆ _enthalpy_varnum

const unsigned int PorousFlowFluidStateSingleComponent::_enthalpy_varnum
protected

Moose variable number of the enthalpy.

Definition at line 62 of file PorousFlowFluidStateSingleComponent.h.

◆ _fluid_density

MaterialProperty<std::vector<Real> >& PorousFlowFluidStateSingleComponent::_fluid_density
protected

Fluid density of each phase.

Definition at line 88 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties(), initQpStatefulProperties(), and setMaterialVectorSize().

◆ _fluid_enthalpy

MaterialProperty<std::vector<Real> >& PorousFlowFluidStateSingleComponent::_fluid_enthalpy
protected

Enthalpy of each phase.

Definition at line 96 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties(), initQpStatefulProperties(), and setMaterialVectorSize().

◆ _fluid_internal_energy

MaterialProperty<std::vector<Real> >& PorousFlowFluidStateSingleComponent::_fluid_internal_energy
protected

Internal energy of each phase.

Definition at line 100 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties(), initQpStatefulProperties(), and setMaterialVectorSize().

◆ _fluid_viscosity

MaterialProperty<std::vector<Real> >& PorousFlowFluidStateSingleComponent::_fluid_viscosity
protected

Viscosity of each phase.

Definition at line 92 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties(), initQpStatefulProperties(), and setMaterialVectorSize().

◆ _fs

const PorousFlowFluidStateSingleComponentBase& PorousFlowFluidStateSingleComponent::_fs
protected

◆ _fsp

std::vector<FluidStateProperties> PorousFlowFluidStateSingleComponent::_fsp
protected

◆ _gas_phase_number

const unsigned int PorousFlowFluidStateSingleComponent::_gas_phase_number
protected

Phase number of the gas phase.

Definition at line 70 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties().

◆ _grad_mass_frac_qp

MaterialProperty<std::vector<std::vector<RealGradient> > >* PorousFlowFluidStateSingleComponent::_grad_mass_frac_qp
protected

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

Definition at line 84 of file PorousFlowFluidStateSingleComponent.h.

◆ _grad_temperature_qp

MaterialProperty<RealGradient>* PorousFlowFluidStateSingleComponent::_grad_temperature_qp
protected

Gradient of temperature (only defined at the qps)

Definition at line 74 of file PorousFlowFluidStateSingleComponent.h.

◆ _gradh_qp

const VariableGradient& PorousFlowFluidStateSingleComponent::_gradh_qp
protected

Gradient of enthalpy (only defined at the qps)

Definition at line 60 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties().

◆ _gradp_qp

MaterialProperty<std::vector<RealGradient> >* const PorousFlowVariableBase::_gradp_qp
protectedinherited

Grad(p) at the quadpoints.

Definition at line 50 of file PorousFlowVariableBase.h.

◆ _grads_qp

MaterialProperty<std::vector<RealGradient> >* const PorousFlowVariableBase::_grads_qp
protectedinherited

Grad(s) at the quadpoints.

Definition at line 65 of file PorousFlowVariableBase.h.

◆ _hidx

const unsigned int PorousFlowFluidStateSingleComponent::_hidx
protected

Index of derivative wrt enthalpy.

Definition at line 116 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties().

◆ _hvar

const unsigned int PorousFlowFluidStateSingleComponent::_hvar
protected

PorousFlow variable number of the enthalpy.

Definition at line 64 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties().

◆ _is_initqp

bool PorousFlowFluidStateSingleComponent::_is_initqp
protected

Flag to indicate whether to calculate stateful properties.

Definition at line 106 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties(), initQpStatefulProperties(), and setMaterialVectorSize().

◆ _liquid_gradp_qp

const VariableGradient& PorousFlowFluidStateSingleComponent::_liquid_gradp_qp
protected

Gradient of porepressure (only defined at the qps)

Definition at line 52 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties().

◆ _liquid_porepressure

const VariableValue& PorousFlowFluidStateSingleComponent::_liquid_porepressure
protected

Porepressure.

Definition at line 50 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties(), and thermophysicalProperties().

◆ _liquid_porepressure_varnum

const unsigned int PorousFlowFluidStateSingleComponent::_liquid_porepressure_varnum
protected

Moose variable number of the porepressure.

Definition at line 54 of file PorousFlowFluidStateSingleComponent.h.

◆ _mass_frac

MaterialProperty<std::vector<std::vector<Real> > >& PorousFlowFluidStateSingleComponent::_mass_frac
protected

Mass fraction matrix.

Definition at line 82 of file PorousFlowFluidStateSingleComponent.h.

Referenced by setMaterialVectorSize().

◆ _num_components

const unsigned int PorousFlowVariableBase::_num_components
protectedinherited

◆ _num_pf_vars

const unsigned int PorousFlowVariableBase::_num_pf_vars
protectedinherited

◆ _num_phases

const unsigned int PorousFlowVariableBase::_num_phases
protectedinherited

◆ _pc

const PorousFlowCapillaryPressure& PorousFlowFluidStateSingleComponent::_pc
protected

Capillary pressure UserObject.

Definition at line 112 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties().

◆ _phase_state

FluidStatePhaseEnum PorousFlowFluidStateSingleComponent::_phase_state
protected

FluidStatePhaseEnum.

Definition at line 110 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties(), and thermophysicalProperties().

◆ _pidx

const unsigned int PorousFlowFluidStateSingleComponent::_pidx
protected

Index of derivative wrt pressure.

Definition at line 114 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties().

◆ _porepressure

MaterialProperty<std::vector<Real> >& PorousFlowVariableBase::_porepressure
protectedinherited

◆ _pvar

const unsigned int PorousFlowFluidStateSingleComponent::_pvar
protected

PorousFlow variable number of the porepressure.

Definition at line 56 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties().

◆ _saturation

MaterialProperty<std::vector<Real> >& PorousFlowVariableBase::_saturation
protectedinherited

◆ _T_c2k

const Real PorousFlowFluidStateSingleComponent::_T_c2k
protected

Conversion from degrees Celsius to degrees Kelvin.

Definition at line 104 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties(), and initQpStatefulProperties().

◆ _temperature

MaterialProperty<Real>& PorousFlowFluidStateSingleComponent::_temperature
protected

Temperature.

Definition at line 72 of file PorousFlowFluidStateSingleComponent.h.

Referenced by computeQpProperties(), and initQpStatefulProperties().


The documentation for this class was generated from the following files:
PorousFlowFluidStateSingleComponent::_phase_state
FluidStatePhaseEnum _phase_state
FluidStatePhaseEnum.
Definition: PorousFlowFluidStateSingleComponent.h:110
PorousFlowVariableBase::computeQpProperties
virtual void computeQpProperties() override
Definition: PorousFlowVariableBase.C:74
PorousFlowFluidStateSingleComponent::_pidx
const unsigned int _pidx
Index of derivative wrt pressure.
Definition: PorousFlowFluidStateSingleComponent.h:114
PorousFlowFluidStateSingleComponent::_fluid_viscosity
MaterialProperty< std::vector< Real > > & _fluid_viscosity
Viscosity of each phase.
Definition: PorousFlowFluidStateSingleComponent.h:92
PorousFlowFluidStateSingleComponent::_dfluid_viscosity_dvar
MaterialProperty< std::vector< std::vector< Real > > > & _dfluid_viscosity_dvar
Derivative of the fluid viscosity for each phase wrt PorousFlow variables.
Definition: PorousFlowFluidStateSingleComponent.h:94
PorousFlowFluidStateSingleComponent::_T_c2k
const Real _T_c2k
Conversion from degrees Celsius to degrees Kelvin.
Definition: PorousFlowFluidStateSingleComponent.h:104
PorousFlowVariableBase::_num_components
const unsigned int _num_components
Number of components.
Definition: PorousFlowVariableBase.h:38
PorousFlowFluidStateSingleComponent::_dgrad_temperature_dgradv
MaterialProperty< std::vector< Real > > *const _dgrad_temperature_dgradv
d(grad temperature)/d(grad PorousFlow variable) at the quadpoints
Definition: PorousFlowFluidStateSingleComponent.h:78
PorousFlowFluidStateSingleComponent::_liquid_porepressure_varnum
const unsigned int _liquid_porepressure_varnum
Moose variable number of the porepressure.
Definition: PorousFlowFluidStateSingleComponent.h:54
libMesh::RealGradient
VectorValue< Real > RealGradient
Definition: GrainForceAndTorqueInterface.h:17
PorousFlowFluidStateSingleComponent::_fsp
std::vector< FluidStateProperties > _fsp
FluidStateProperties data structure.
Definition: PorousFlowFluidStateSingleComponent.h:108
PorousFlowCapillaryPressure
Base class for capillary pressure for multiphase flow in porous media.
Definition: PorousFlowCapillaryPressure.h:39
PorousFlowFluidStateSingleComponent::thermophysicalProperties
virtual void thermophysicalProperties()
Calculates all required thermophysical properties and derivatives for each phase and fluid component.
Definition: PorousFlowFluidStateSingleComponent.C:136
PorousFlowFluidStateSingleComponent::_dmass_frac_dvar
MaterialProperty< std::vector< std::vector< std::vector< Real > > > > & _dmass_frac_dvar
Derivative of the mass fraction matrix with respect to the Porous Flow variables.
Definition: PorousFlowFluidStateSingleComponent.h:86
PorousFlowFluidStateSingleComponent::_pvar
const unsigned int _pvar
PorousFlow variable number of the porepressure.
Definition: PorousFlowFluidStateSingleComponent.h:56
PorousFlowFluidStateSingleComponent::_is_initqp
bool _is_initqp
Flag to indicate whether to calculate stateful properties.
Definition: PorousFlowFluidStateSingleComponent.h:106
PorousFlowFluidStateSingleComponent::_fluid_internal_energy
MaterialProperty< std::vector< Real > > & _fluid_internal_energy
Internal energy of each phase.
Definition: PorousFlowFluidStateSingleComponent.h:100
PorousFlowFluidStateSingleComponent::_temperature
MaterialProperty< Real > & _temperature
Temperature.
Definition: PorousFlowFluidStateSingleComponent.h:72
PorousFlowFluidStateSingleComponentBase
Base class for miscible multiphase flow classes with a single fluid component using a pressure and en...
Definition: PorousFlowFluidStateSingleComponentBase.h:23
PorousFlowVariableBase::_num_phases
const unsigned int _num_phases
Number of phases.
Definition: PorousFlowVariableBase.h:35
PorousFlowVariableBase::_dgradp_qp_dgradv
MaterialProperty< std::vector< std::vector< Real > > > *const _dgradp_qp_dgradv
d(grad porepressure)/d(grad PorousFlow variable) at the quadpoints
Definition: PorousFlowVariableBase.h:53
PorousFlowFluidStateBase::numPhases
unsigned int numPhases() const
The maximum number of phases in this model.
Definition: PorousFlowFluidStateBase.h:68
PorousFlowFluidStateSingleComponent::_dfluid_enthalpy_dvar
MaterialProperty< std::vector< std::vector< Real > > > & _dfluid_enthalpy_dvar
Derivative of the fluid enthalpy for each phase wrt PorousFlow variables.
Definition: PorousFlowFluidStateSingleComponent.h:98
PorousFlowFluidStateSingleComponent::_liquid_gradp_qp
const VariableGradient & _liquid_gradp_qp
Gradient of porepressure (only defined at the qps)
Definition: PorousFlowFluidStateSingleComponent.h:52
PorousFlowFluidStateSingleComponent::_pc
const PorousFlowCapillaryPressure & _pc
Capillary pressure UserObject.
Definition: PorousFlowFluidStateSingleComponent.h:112
PorousFlowFluidStateSingleComponent::_dgrad_temperature_dv
MaterialProperty< std::vector< RealGradient > > *const _dgrad_temperature_dv
d(grad temperature)/d(PorousFlow variable) at the quadpoints
Definition: PorousFlowFluidStateSingleComponent.h:80
PorousFlowFluidStateSingleComponent::_enthalpy_varnum
const unsigned int _enthalpy_varnum
Moose variable number of the enthalpy.
Definition: PorousFlowFluidStateSingleComponent.h:62
PorousFlowFluidStateSingleComponentBase::thermophysicalProperties
virtual void thermophysicalProperties(Real pressure, Real enthalpy, unsigned int qp, FluidStatePhaseEnum &phase_state, std::vector< FluidStateProperties > &fsp) const =0
Determines the complete thermophysical state of the system for a given set of primary variables.
PorousFlowFluidStateSingleComponent::_fs
const PorousFlowFluidStateSingleComponentBase & _fs
FluidState UserObject.
Definition: PorousFlowFluidStateSingleComponent.h:66
FluidStateProperties
AD data structure to pass calculated thermophysical properties.
Definition: PorousFlowFluidStateBase.h:26
PorousFlowFluidStateSingleComponent::_dtemperature_dvar
MaterialProperty< std::vector< Real > > & _dtemperature_dvar
Derivative of temperature wrt PorousFlow variables.
Definition: PorousFlowFluidStateSingleComponent.h:76
PorousFlowVariableBase::_porepressure
MaterialProperty< std::vector< Real > > & _porepressure
Computed nodal or quadpoint values of porepressure of the phases.
Definition: PorousFlowVariableBase.h:44
PorousFlowFluidStateSingleComponent::_aqueous_phase_number
const unsigned int _aqueous_phase_number
Phase number of the aqueous phase.
Definition: PorousFlowFluidStateSingleComponent.h:68
PorousFlowVariableBase::_dporepressure_dvar
MaterialProperty< std::vector< std::vector< Real > > > & _dporepressure_dvar
d(porepressure)/d(PorousFlow variable)
Definition: PorousFlowVariableBase.h:47
PorousFlowVariableBase::_num_pf_vars
const unsigned int _num_pf_vars
Number of PorousFlow variables.
Definition: PorousFlowVariableBase.h:41
name
const std::string name
Definition: Setup.h:21
PorousFlowFluidStateSingleComponent::_dfluid_density_dvar
MaterialProperty< std::vector< std::vector< Real > > > & _dfluid_density_dvar
Derivative of the fluid density for each phase wrt PorousFlow variables.
Definition: PorousFlowFluidStateSingleComponent.h:90
PorousFlowFluidStateSingleComponent::_enthalpy
const VariableValue & _enthalpy
Enthalpy.
Definition: PorousFlowFluidStateSingleComponent.h:58
PorousFlowFluidStateSingleComponent::_gas_phase_number
const unsigned int _gas_phase_number
Phase number of the gas phase.
Definition: PorousFlowFluidStateSingleComponent.h:70
PorousFlowCapillaryPressure::dCapillaryPressure
virtual Real dCapillaryPressure(Real saturation, unsigned qp=0) const
Derivative of capillary pressure wrt true saturation.
Definition: PorousFlowCapillaryPressure.C:73
PorousFlowFluidStateSingleComponent::_hidx
const unsigned int _hidx
Index of derivative wrt enthalpy.
Definition: PorousFlowFluidStateSingleComponent.h:116
PorousFlowFluidStateBase::clearFluidStateProperties
void clearFluidStateProperties(std::vector< FluidStateProperties > &fsp) const
Clears the contents of the FluidStateProperties data structure.
Definition: PorousFlowFluidStateBase.C:39
PorousFlowFluidStateSingleComponent::_dfluid_internal_energy_dvar
MaterialProperty< std::vector< std::vector< Real > > > & _dfluid_internal_energy_dvar
Derivative of the fluid internal energy for each phase wrt PorousFlow variables.
Definition: PorousFlowFluidStateSingleComponent.h:102
PorousFlowFluidStateSingleComponent::_fluid_enthalpy
MaterialProperty< std::vector< Real > > & _fluid_enthalpy
Enthalpy of each phase.
Definition: PorousFlowFluidStateSingleComponent.h:96
PorousFlowFluidStateSingleComponent::_gradh_qp
const VariableGradient & _gradh_qp
Gradient of enthalpy (only defined at the qps)
Definition: PorousFlowFluidStateSingleComponent.h:60
PorousFlowFluidStateSingleComponent::_hvar
const unsigned int _hvar
PorousFlow variable number of the enthalpy.
Definition: PorousFlowFluidStateSingleComponent.h:64
PorousFlowCapillaryPressure::d2CapillaryPressure
virtual Real d2CapillaryPressure(Real saturation, unsigned qp=0) const
Second derivative of capillary pressure wrt true saturation.
Definition: PorousFlowCapillaryPressure.C:82
PorousFlowFluidStateSingleComponent::setMaterialVectorSize
void setMaterialVectorSize() const
Size material property vectors and initialise with zeros.
Definition: PorousFlowFluidStateSingleComponent.C:341
PorousFlowFluidStateSingleComponent::_mass_frac
MaterialProperty< std::vector< std::vector< Real > > > & _mass_frac
Mass fraction matrix.
Definition: PorousFlowFluidStateSingleComponent.h:82
PorousFlowVariableBase::PorousFlowVariableBase
PorousFlowVariableBase(const InputParameters &parameters)
Definition: PorousFlowVariableBase.C:23
PorousFlowFluidStateSingleComponent::_grad_temperature_qp
MaterialProperty< RealGradient > * _grad_temperature_qp
Gradient of temperature (only defined at the qps)
Definition: PorousFlowFluidStateSingleComponent.h:74
PorousFlowFluidStateSingleComponent::_liquid_porepressure
const VariableValue & _liquid_porepressure
Porepressure.
Definition: PorousFlowFluidStateSingleComponent.h:50
PorousFlowVariableBase::initQpStatefulProperties
virtual void initQpStatefulProperties() override
Definition: PorousFlowVariableBase.C:66
PorousFlowVariableBase::_saturation
MaterialProperty< std::vector< Real > > & _saturation
Computed nodal or qp saturation of the phases.
Definition: PorousFlowVariableBase.h:59
PorousFlowVariableBase::_dsaturation_dvar
MaterialProperty< std::vector< std::vector< Real > > > & _dsaturation_dvar
d(saturation)/d(PorousFlow variable)
Definition: PorousFlowVariableBase.h:62
PorousFlowFluidStateSingleComponent::_grad_mass_frac_qp
MaterialProperty< std::vector< std::vector< RealGradient > > > * _grad_mass_frac_qp
Gradient of the mass fraction matrix (only defined at the qps)
Definition: PorousFlowFluidStateSingleComponent.h:84
PorousFlowFluidStateSingleComponent::_fluid_density
MaterialProperty< std::vector< Real > > & _fluid_density
Fluid density of each phase.
Definition: PorousFlowFluidStateSingleComponent.h:88