- porepressureThe name of the porepressure variable
C++ Type:VariableName
Controllable:No
Description:The name of the porepressure variable
PorousFlowFullySaturated
Adds Kernels and fluid-property Materials necessary to simulate a single-phase fully-saturated flow problem. Numerical stabilization options for the fluid and heat flow are: no upwinding, full-upwinding or KT stabilization. No Kernels for diffusion and dispersion of fluid components are added. To run a simulation you will also need to provide various other Materials for each mesh block, depending on your simulation type, viz: permeability, porosity, elasticity tensor, strain calculator, stress calculator, matrix internal energy, thermal conductivity, diffusivity
For a worked example, see tutorial page 6.
Dictator
This Action
automatically adds an appropriately-parameterized PorousFlowDictator with name given by the dictator_name
input.
Fluid equations
Order of unknowns
The unknowns for the fluid equations are Here is the number of fluid components, is a fluid component, and is the porepressure. For example:
For single-component fluids () such as pure water, is the only variable
For brine, whose fluid components are NaCl and HO, the variables are
For an acidic brine, whose basis species are Na, Cl, H and HO (with ), the variables are
Note that this Action
associates the final mass-fraction with the porepressure variable. Ordering your mass-fraction variables so that the final variable has the greatest mass-fraction will usually improve convergence speed because the Jacobian is more diagonally dominant. For instance, for the brine situation, it is usually much more computationally efficient to choose rather than even though the final answers will be identical.
Since the final mass-fraction, , is associated with the porepressure variable, BCs
and other objects associated with the porepressure variable will be "acting on" the final mass-fraction. See tutorial page 6 for examples.
Fluid-flow equations
This Action
simulates the fluid equations (1) In this equation, the fluid density, , appears in parenthases, because the user has the option of including it or not using the multiply_by_density
flag.
Think carefully before you use multiply_by_density = false
. When not multiply by density, in many situations the sum of Eq. (1) will imply a Laplace equation such , that is a steady-state distribution for the porepressure. In addition, care must be taken when using other parts of PorousFlow, for instance, the PorousFlowMass
Postprocessor is coded to record fluid mass, not fluid volume. New users should set multiply_by_density = true
to avoid confusion, even at the expense of extra computational cost.
Kernels added
To represent the term (which only appears in transient
simulations) the Action
adds a PorousFlowMasstimeDerivative Kernel for each fluid component. These Kernels lump the fluid-component mass to the nodes to ensure superior numerical stabilization.
To represent the term (which only appears in transient
simulations with mechanical coupling_type
) the Action
adds a PorousFlowMassVolumetricExpansion Kernel for each fluid component. These Kernels lump the fluid-component mass to the nodes to ensure superior numerical stabilization.
To represent the the Action
adds Kernels and possibly UserObjects depending on the stabilization
parameter:
In the case of no-upwinding, a PorousFlowFullySaturatedDarcyFlow Kernel for each fluid component.
For full upwinding, a PorousFlowFullySaturatedAdvectiveFlux Kernel for each fluid component.
For KT stabilization, a PorousFlowFluxLimitedTVDAdvection Kernel and a PorousFlowAdvectiveFluxCalculatorSaturatedMultiComponent UserObject for each fluid component (or a PorousFlowAdvectiveFluxCalculatorSaturated in the case of just one fluid component).
Advanced additions: diffusion and dispersion, radioactive decay, chemistry
The Action
does not add Kernels to model the following effects
Dispersion: if dispersion is important, the user must manually add a set of PorousFlowDispersiveFlux Kernels.
Radioactive decay: if this is important, the user must manually add some PorousFlowMassRadioactiveDecay Kernels.
Chemistry: if chemistry is of importance, then appropriate Kernels must be manually added, such as PorousFlowPreDis, PorousFlowDesorpedMassTimeDerivative and PorousFlowDesorpedMassVolumetricExpansion. A simple example is discussed in the tutorial 07 page. Alternatively, MOOSE's Geochemistry module may be employed.
Heat equation
The heat equation is associated with the temperature variable, and is only included if the coupling_type
includes "Thermo". (2) Note that multiplication by fluid density occurs, irrespective of the multiply_by_density
flag.
To represent the term (which only appears in transient
simulations) the Action
adds a PorousFlowEnergyTimeDerivative Kernel. This Kernel lumps the heat energy-density to the nodes.
To represent the term (which only appears in transient
simulations with coupling_type = ThermoHydroMechanical
), the Action
adds a PorousFlowHeatVolumetricExpansion Kernel.
To represent the term , the Action
adds a PorousFlowHeatConduction Kernel
To represent the term , the Action
adds the following Kernels and possibly UserObjects, depending on the stabilization
parameters:
In the case of no-upwinding, a PorousFlowFullySaturatedHeatAdvection Kernel with
multiply_by_density = true
.For full upwinding, a PorousFlowFullySaturatedUpwindHeatAdvection Kernel.
For KT stabilization, a PorousFlowFluxLimitedTVDAdvection Kernel and a PorousFlowAdvectiveFluxCalculatorSaturatedHeat UserObject with
multiply_by_density = true
.
Solid-mechanics equations
The static conservation of momentum equation is associated with the displacement variables, and is only included if the coupling_type
includes "Mechanical": (3) These equations are represented by a set of StressDivergenceTensors
Kernels (or StressDivergenceRZTensors
Kernels in the case of RZ coordinates), including thermal-eigenstrain coupling if coupling_type = HydroThermoMechanical
, and Gravity
Kernels. These Kernels are added automatically and are part of the Tensor-Mechanics module. The term is represented by a PorousFlowEffectiveStressCoupling Kernel.
It is assumed that the effective stress not the total stress enters into the consitutive law (as above), and any plasticity, and any insitu stresses, and almost everywhere else. One exception is specifying Neumann boundary conditions for the displacements where the total stresses are being specified, as can be seen from Eq. (3). Therefore, MOOSE uses effective stress, and not total stress, internally. If one needs to input or output total stress, one must subtract from MOOSE's stress.
Sources and boundary conditions
Source terms may be added to the above equations by the user by including appropriate Kernels
or DiracKernels
(such as those listed here and described here). Boundary conditions can also be suppled. Sources and boundary conditions are not automatically included by the Action
.
Materials added
If required by the Kernels, the following Materials are added. For each case, derivatives of the Material Properties are also created: these are with respect to mass fractions and porepressure, and, in the case of appropriate coupling_type
, temperature and displacements.
PorousFlow1PhaseFullySaturated which creates Material
versions of the porepressure and its spatial derivatives, saturation (which is always 1.0 in this case) and its spatial derivatives.
PorousFlowTemperature which creates Material
versions of the temperature and its spatial derivatives.
PorousFlowMassFraction which creates Material
versions of the mass-fraction variables and their spatial derivatives
PorousFlowSingleComponentFluid or PorousFlowBrine or, in the case of custom_fluid_properties
, a user-supplied Material to compute the fluid properties. These properties are the fluid density and viscosity, as well as internal energy and enthalpy if required.
Using non-standard pressure_unit
or time_unit
is only possible with PorousFlowSingleComponentFluid
, not PorousFlowBrine
. You need to read the essay on these non-standard unit systems here
PorousFlowEffectiveFluidPressure which creates an effective fluid pressure, equal to the porepressure, for use in mechanically-coupled models.
PorousFlowVolumetricStrain which creates a volumetric strain rate, for use in mechanically-coupled models.
PorousFlowNearestQp which records the nearest quadpoint to each node in each element.
Even though there is only one fluid phase in this fully saturated action, some objects may require a relative permeability material to work. Examples include PorousFlowDarcyVelocityComponent AuxKernels
which requires relative permeability, or PorousFlowPeacemanBorehole wells. By default, this action adds a constant relative permeability of one, so that the fluid is perfectly mobile.
Materials not added
Various important Materials
are not added by this Action, so must be added by the user in the [Materials]
block. The reason these are not added by default is that they are usually subdomain-dependent.
One of the PorousFlowPorosity options
One of the PorousFlowPermeability options
Thermal conductivity, such as PorousFlowThermalConductivityIdeal
A definition of the elasticity tensor (eg ComputeIsotropicElasticityTensor) a strain calculator (eg ComputeSmallStrain) a thermal expansion eigenstrain calculator, (eg ComputeThermalExpansionEigenstrain) and a stress calculator (eg ComputeLinearElasticStress)
Materials associated with chemistry, such as PorousFlowAqueousPreDisChemistry and PorousFlowAqueousPreDisMineral: see tutorial 07 for a worked example.
AuxVariables
If
add_darcy_aux = true
then theAction
addsconstant monomial
auxillary variables with the namesdarcy_vel_x
,darcy_vel_y
anddarcy_vel_z
that record the Darcy velocity by using PorousFlowDarcyVelocityComponent AuxKernels.If
add_stress_aux = true
and thecoupling_type
includes "Mechanical", then theAction
addsconstant monomial
auxillary variables with the namestress_xx
,stress_xy
,stress_xz
,stress_yx
,stress_yy
,stress_yz
,stress_zx
,stress_zy
andstress_zz
that record the effective stress using aRankTwoAux
. As mentioned above, this is the effective stress: if you require total stress you need to subtract , as in
Input Parameters
- active__all__ If specified only the blocks named will be visited and made active
Default:__all__
C++ Type:std::vector<std::string>
Controllable:No
Description:If specified only the blocks named will be visited and made active
- add_darcy_auxTrueAdd AuxVariables that record Darcy velocity
Default:True
C++ Type:bool
Controllable:No
Description:Add AuxVariables that record Darcy velocity
- add_stress_auxTrueAdd AuxVariables that record effective stress
Default:True
C++ Type:bool
Controllable:No
Description:Add AuxVariables that record effective stress
- base_nameThe base_name used in the TensorMechanics strain calculator. This is only relevant for mechanically-coupled models.
C++ Type:std::string
Controllable:No
Description:The base_name used in the TensorMechanics strain calculator. This is only relevant for mechanically-coupled models.
- biot_coefficient1The Biot coefficient (relevant only for mechanically-coupled simulations)
Default:1
C++ Type:double
Controllable:No
Description:The Biot coefficient (relevant only for mechanically-coupled simulations)
- coupling_typeHydroThe type of simulation. For simulations involving Mechanical deformations, you will need to supply the correct Biot coefficient. For simulations involving Thermal flows, you will need an associated ConstantThermalExpansionCoefficient Material
Default:Hydro
C++ Type:MooseEnum
Controllable:No
Description:The type of simulation. For simulations involving Mechanical deformations, you will need to supply the correct Biot coefficient. For simulations involving Thermal flows, you will need an associated ConstantThermalExpansionCoefficient Material
- dictator_namedictatorThe name of the dictator user object that is created by this Action
Default:dictator
C++ Type:std::string
Controllable:No
Description:The name of the dictator user object that is created by this Action
- displacementsThe name of the displacement variables (relevant only for mechanically-coupled simulations)
C++ Type:std::vector<VariableName>
Controllable:No
Description:The name of the displacement variables (relevant only for mechanically-coupled simulations)
- eigenstrain_namesList of all eigenstrain models used in mechanics calculations. Typically the eigenstrain_name used in ComputeThermalExpansionEigenstrain. Only needed for thermally-coupled simulations with thermal expansion.
C++ Type:std::vector<MaterialPropertyName>
Controllable:No
Description:List of all eigenstrain models used in mechanics calculations. Typically the eigenstrain_name used in ComputeThermalExpansionEigenstrain. Only needed for thermally-coupled simulations with thermal expansion.
- fluid_properties_typePorousFlowSingleComponentFluidType of fluid properties to use. For 'PorousFlowSingleComponentFluid' you must provide a fp UserObject. For 'PorousFlowBrine' you must supply a nacl_name. For 'Custom' your input file must include a Material that provides fluid properties such as density, viscosity, enthalpy and internal energy
Default:PorousFlowSingleComponentFluid
C++ Type:MooseEnum
Controllable:No
Description:Type of fluid properties to use. For 'PorousFlowSingleComponentFluid' you must provide a fp UserObject. For 'PorousFlowBrine' you must supply a nacl_name. For 'Custom' your input file must include a Material that provides fluid properties such as density, viscosity, enthalpy and internal energy
- flux_limiter_typeVanLeerType of flux limiter to use if stabilization=KT. 'None' means that no antidiffusion will be added in the Kuzmin-Turek scheme
Default:VanLeer
C++ Type:MooseEnum
Controllable:No
Description:Type of flux limiter to use if stabilization=KT. 'None' means that no antidiffusion will be added in the Kuzmin-Turek scheme
- fpThe name of the user object for fluid properties. Only needed if fluid_properties_type = PorousFlowSingleComponentFluid
C++ Type:UserObjectName
Controllable:No
Description:The name of the user object for fluid properties. Only needed if fluid_properties_type = PorousFlowSingleComponentFluid
- gravity0 0 -10Gravitational acceleration vector downwards (m/s^2)
Default:0 0 -10
C++ Type:libMesh::VectorValue<double>
Controllable:No
Description:Gravitational acceleration vector downwards (m/s^2)
- inactiveIf specified blocks matching these identifiers will be skipped.
C++ Type:std::vector<std::string>
Controllable:No
Description:If specified blocks matching these identifiers will be skipped.
- mass_fraction_varsList of variables that represent the mass fractions. With only one fluid component, this may be left empty. With N fluid components, the format is 'f_0 f_1 f_2 ... f_(N-1)'. That is, the N^th component need not be specified because f_N = 1 - (f_0 + f_1 + ... + f_(N-1)). It is best numerically to choose the N-1 mass fraction variables so that they represent the fluid components with small concentrations. This Action will associated the i^th mass fraction variable to the equation for the i^th fluid component, and the pressure variable to the N^th fluid component.
C++ Type:std::vector<VariableName>
Controllable:No
Description:List of variables that represent the mass fractions. With only one fluid component, this may be left empty. With N fluid components, the format is 'f_0 f_1 f_2 ... f_(N-1)'. That is, the N^th component need not be specified because f_N = 1 - (f_0 + f_1 + ... + f_(N-1)). It is best numerically to choose the N-1 mass fraction variables so that they represent the fluid components with small concentrations. This Action will associated the i^th mass fraction variable to the equation for the i^th fluid component, and the pressure variable to the N^th fluid component.
- multiply_by_densityTrueIf true, then the Kernels for fluid flow are multiplied by the fluid density. If false, this multiplication is not performed, which means the problem becomes more linear, but care must be taken when using other PorousFlow objects, since MOOSE will be computing volume fluxes, not mass fluxes.
Default:True
C++ Type:bool
Controllable:No
Description:If true, then the Kernels for fluid flow are multiplied by the fluid density. If false, this multiplication is not performed, which means the problem becomes more linear, but care must be taken when using other PorousFlow objects, since MOOSE will be computing volume fluxes, not mass fluxes.
- nacl_nameName of the NaCl variable. Only required if fluid_properties_type = PorousFlowBrine
C++ Type:VariableName
Controllable:No
Description:Name of the NaCl variable. Only required if fluid_properties_type = PorousFlowBrine
- number_aqueous_equilibrium0The number of secondary species in the aqueous-equilibrium reaction system. (Leave as zero if the simulation does not involve chemistry)
Default:0
C++ Type:unsigned int
Controllable:No
Description:The number of secondary species in the aqueous-equilibrium reaction system. (Leave as zero if the simulation does not involve chemistry)
- number_aqueous_kinetic0The number of secondary species in the aqueous-kinetic reaction system involved in precipitation and dissolution. (Leave as zero if the simulation does not involve chemistry)
Default:0
C++ Type:unsigned int
Controllable:No
Description:The number of secondary species in the aqueous-kinetic reaction system involved in precipitation and dissolution. (Leave as zero if the simulation does not involve chemistry)
- pressure_unitPaThe unit of the pressure variable used everywhere in the input file except for in the FluidProperties-module objects. This can be set to the non-default value only for fluid_properties_type = PorousFlowSingleComponentFluid
Default:Pa
C++ Type:MooseEnum
Controllable:No
Description:The unit of the pressure variable used everywhere in the input file except for in the FluidProperties-module objects. This can be set to the non-default value only for fluid_properties_type = PorousFlowSingleComponentFluid
- save_component_rate_inList of AuxVariables into which the rate-of-change of each fluid component at each node will be saved. There must be exactly N of these to match the N fluid components. The result will be measured in kg/s, where the kg is the mass of the fluid component at the node (or m^3/s if multiply_by_density=false). Note that this saves the result from the MassTimeDerivative Kernels, but NOT from the MassVolumetricExpansion Kernels.
C++ Type:std::vector<AuxVariableName>
Controllable:No
Description:List of AuxVariables into which the rate-of-change of each fluid component at each node will be saved. There must be exactly N of these to match the N fluid components. The result will be measured in kg/s, where the kg is the mass of the fluid component at the node (or m^3/s if multiply_by_density=false). Note that this saves the result from the MassTimeDerivative Kernels, but NOT from the MassVolumetricExpansion Kernels.
- stabilizationFullNumerical stabilization used. 'Full' means full upwinding. 'KT' means FEM-TVD stabilization of Kuzmin-Turek
Default:Full
C++ Type:MooseEnum
Controllable:No
Description:Numerical stabilization used. 'Full' means full upwinding. 'KT' means FEM-TVD stabilization of Kuzmin-Turek
- strain_at_nearest_qpFalseOnly relevant for models in which porosity depends on strain. If true, then when calculating nodal porosity that depends on strain, the strain at the nearest quadpoint will be used. This adds a small extra computational burden, and is only necessary for simulations involving: (1) elements that are not linear lagrange or (2) certain PorousFlow Dirac Kernels (as specified in their documentation). If you set this to true, you will also want to set the same parameter to true for related Kernels and Materials (which is probably easiest to do in the GlobalParams block)
Default:False
C++ Type:bool
Controllable:No
Description:Only relevant for models in which porosity depends on strain. If true, then when calculating nodal porosity that depends on strain, the strain at the nearest quadpoint will be used. This adds a small extra computational burden, and is only necessary for simulations involving: (1) elements that are not linear lagrange or (2) certain PorousFlow Dirac Kernels (as specified in their documentation). If you set this to true, you will also want to set the same parameter to true for related Kernels and Materials (which is probably easiest to do in the GlobalParams block)
- temperature293.0For isothermal simulations, this is the temperature at which fluid properties (and stress-free strains) are evaluated at. Otherwise, this is the name of the temperature variable. Units = Kelvin
Default:293.0
C++ Type:std::vector<VariableName>
Controllable:No
Description:For isothermal simulations, this is the temperature at which fluid properties (and stress-free strains) are evaluated at. Otherwise, this is the name of the temperature variable. Units = Kelvin
- temperature_unitKelvinThe unit of the temperature variable
Default:Kelvin
C++ Type:MooseEnum
Controllable:No
Description:The unit of the temperature variable
- time_unitsecondsThe unit of time used everywhere in the input file except for in the FluidProperties-module objects. This can be set to the non-default value only for fluid_properties_type = PorousFlowSingleComponentFluid
Default:seconds
C++ Type:MooseEnum
Controllable:No
Description:The unit of time used everywhere in the input file except for in the FluidProperties-module objects. This can be set to the non-default value only for fluid_properties_type = PorousFlowSingleComponentFluid
- use_displaced_meshFalseUse displaced mesh computations in mechanical kernels
Default:False
C++ Type:bool
Controllable:No
Description:Use displaced mesh computations in mechanical kernels