Mass computation and conservation
The total fluid mass of species within a volume is (1) It must be checked that MOOSE calculates this correctly, using the PorousFlowFluidMass postprocessor, in order that mass-balances be correct, and also because this quantity is used in a number of other tests.
Single-phase, single-component fluid
A 1D model with , and with three elements of size 1 is created with the following properties:
Table 1: Parameter values in the single-phase, single-component test
Parameter | Value |
---|---|
Constant fluid bulk modulus | Pa |
Fluid density at zero pressure | kg.m |
Van Genuchten | 0.5 |
Van Genuchten | Pa |
Porosity | 0.1 |
The porepressure is set at .
# checking that the mass postprocessor correctly calculates the mass
# 1phase, 1component, constant porosity
[Mesh<<<{"href": "../../../../syntax/Mesh/index.html"}>>>]
type = GeneratedMesh
dim = 1
nx = 3
xmin = -1
xmax = 1
[]
[GlobalParams<<<{"href": "../../../../syntax/GlobalParams/index.html"}>>>]
PorousFlowDictator = dictator
[]
[Variables<<<{"href": "../../../../syntax/Variables/index.html"}>>>]
[pp]
[]
[]
[ICs<<<{"href": "../../../../syntax/ICs/index.html"}>>>]
[pinit]
type = FunctionIC<<<{"description": "An initial condition that uses a normal function of x, y, z to produce values (and optionally gradients) for a field variable.", "href": "../../../../source/ics/FunctionIC.html"}>>>
function<<<{"description": "The initial condition function."}>>> = x
variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = pp
[]
[]
[Kernels<<<{"href": "../../../../syntax/Kernels/index.html"}>>>]
[mass0]
type = PorousFlowMassTimeDerivative<<<{"description": "Derivative of fluid-component mass with respect to time. Mass lumping to the nodes is used.", "href": "../../../../source/kernels/PorousFlowMassTimeDerivative.html"}>>>
fluid_component<<<{"description": "The index corresponding to the component for this kernel"}>>> = 0
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pp
[]
[]
[UserObjects<<<{"href": "../../../../syntax/UserObjects/index.html"}>>>]
[dictator]
type = PorousFlowDictator<<<{"description": "Holds information on the PorousFlow variable names", "href": "../../../../source/userobjects/PorousFlowDictator.html"}>>>
porous_flow_vars<<<{"description": "List of primary variables that are used in the PorousFlow simulation. Jacobian entries involving derivatives wrt these variables will be computed. In single-phase models you will just have one (eg 'pressure'), in two-phase models you will have two (eg 'p_water p_gas', or 'p_water s_water'), etc."}>>> = 'pp'
number_fluid_phases<<<{"description": "The number of fluid phases in the simulation"}>>> = 1
number_fluid_components<<<{"description": "The number of fluid components in the simulation"}>>> = 1
[]
[pc]
type = PorousFlowCapillaryPressureVG<<<{"description": "van Genuchten capillary pressure", "href": "../../../../source/userobjects/PorousFlowCapillaryPressureVG.html"}>>>
m<<<{"description": "van Genuchten exponent m. Must be between 0 and 1, and optimally should be set to >0.5"}>>> = 0.5
alpha<<<{"description": "van Genuchten parameter alpha. Must be positive"}>>> = 1
[]
[]
[FluidProperties<<<{"href": "../../../../syntax/FluidProperties/index.html"}>>>]
[simple_fluid]
type = SimpleFluidProperties<<<{"description": "Fluid properties for a simple fluid with a constant bulk density", "href": "../../../../source/fluidproperties/SimpleFluidProperties.html"}>>>
bulk_modulus<<<{"description": "Constant bulk modulus (Pa)"}>>> = 1
density0<<<{"description": "Density at zero pressure and zero temperature"}>>> = 1
thermal_expansion<<<{"description": "Constant coefficient of thermal expansion (1/K)"}>>> = 0
[]
[]
[Materials<<<{"href": "../../../../syntax/Materials/index.html"}>>>]
[temperature]
type = PorousFlowTemperature<<<{"description": "Material to provide temperature at the quadpoints or nodes and derivatives of it with respect to the PorousFlow variables", "href": "../../../../source/materials/PorousFlowTemperature.html"}>>>
[]
[ppss]
type = PorousFlow1PhaseP<<<{"description": "This Material is used for the partially saturated single-phase situation where porepressure is the primary variable", "href": "../../../../source/materials/PorousFlow1PhaseP.html"}>>>
porepressure<<<{"description": "Variable that represents the porepressure of the single phase"}>>> = pp
capillary_pressure<<<{"description": "Name of the UserObject defining the capillary pressure"}>>> = pc
[]
[massfrac]
type = PorousFlowMassFraction<<<{"description": "This Material forms a std::vector<std::vector ...> of mass-fractions out of the individual mass fractions", "href": "../../../../source/materials/PorousFlowMassFraction.html"}>>>
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid<<<{"description": "This Material calculates fluid properties at the quadpoints or nodes for a single component fluid", "href": "../../../../source/materials/PorousFlowSingleComponentFluid.html"}>>>
fp<<<{"description": "The name of the user object for fluid properties"}>>> = simple_fluid
phase<<<{"description": "The phase number"}>>> = 0
[]
[porosity]
type = PorousFlowPorosityConst<<<{"description": "This Material calculates the porosity assuming it is constant", "href": "../../../../source/materials/PorousFlowPorosityConst.html"}>>>
porosity<<<{"description": "The porosity (assumed indepenent of porepressure, temperature, strain, etc, for this material). This should be a real number, or a constant monomial variable (not a linear lagrange or other kind of variable)."}>>> = 0.1
[]
[]
[Postprocessors<<<{"href": "../../../../syntax/Postprocessors/index.html"}>>>]
[total_mass]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
[]
[]
[Preconditioning<<<{"href": "../../../../syntax/Preconditioning/index.html"}>>>]
[andy]
type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../../../../source/preconditioners/SingleMatrixPreconditioner.html"}>>>
full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
petsc_options_iname<<<{"description": "Names of PETSc name/value pairs"}>>> = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it'
petsc_options_value<<<{"description": "Values of PETSc name/value pairs (must correspond with \"petsc_options_iname\""}>>> = 'bcgs bjacobi 1 .999 10000'
[]
[]
[Executioner<<<{"href": "../../../../syntax/Executioner/index.html"}>>>]
type = Transient
solve_type = Newton
dt = 1
end_time = 1
[]
[Outputs<<<{"href": "../../../../syntax/Outputs/index.html"}>>>]
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'timestep_end'
file_base<<<{"description": "Common file base name to be utilized with all output objects"}>>> = mass01
csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
[]
(moose/modules/porous_flow/test/tests/mass_conservation/mass01.i)Recall that in PorousFlow, mass is lumped to the nodes. Therefore, the integral above is evaluated at the nodes, and a sum of the results is outputted as the PorousFlowFluidMass postprocessor. Using the properties given above, this yields:
Table 2: Nodal masses in the single-phase, single-component test
Density | Saturation | Nodal mass | ||
---|---|---|---|---|
-1 | -1 | 0.367879441 | 0.707106781 | 0.008671002 |
-0.333333333 | -0.333333333 | 0.716531311 | 0.948683298 | 0.02265871 |
-0.333333333 | -0.333333333 | 0.716531311 | 0.948683298 | 0.02265871 |
0.333333333 | 0.333333333 | 1.395612425 | 1 | 0.046520414 |
0.333333333 | 0.333333333 | 1.395612425 | 1 | 0.046520414 |
1 | 1 | 2.718281828 | 1 | 0.090609394 |
- | - | - | Total | 0.237638643 |
MOOSE also gives the total mass as 0.237638643\,kg.
Single-phase, two-components
This is similar to the previous section, but has two fluid components. The mass fraction is fixed at
# checking that the mass postprocessor correctly calculates the mass
# 1phase, 2component, constant porosity
[Mesh<<<{"href": "../../../../syntax/Mesh/index.html"}>>>]
type = GeneratedMesh
dim = 1
nx = 3
xmin = -1
xmax = 1
[]
[GlobalParams<<<{"href": "../../../../syntax/GlobalParams/index.html"}>>>]
PorousFlowDictator = dictator
[]
[Variables<<<{"href": "../../../../syntax/Variables/index.html"}>>>]
[pp]
[]
[mass_frac_comp0]
[]
[]
[ICs<<<{"href": "../../../../syntax/ICs/index.html"}>>>]
[pinit]
type = FunctionIC<<<{"description": "An initial condition that uses a normal function of x, y, z to produce values (and optionally gradients) for a field variable.", "href": "../../../../source/ics/FunctionIC.html"}>>>
function<<<{"description": "The initial condition function."}>>> = x
variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = pp
[]
[minit]
type = FunctionIC<<<{"description": "An initial condition that uses a normal function of x, y, z to produce values (and optionally gradients) for a field variable.", "href": "../../../../source/ics/FunctionIC.html"}>>>
function<<<{"description": "The initial condition function."}>>> = 'x*x'
variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = mass_frac_comp0
[]
[]
[Kernels<<<{"href": "../../../../syntax/Kernels/index.html"}>>>]
[mass0]
type = PorousFlowMassTimeDerivative<<<{"description": "Derivative of fluid-component mass with respect to time. Mass lumping to the nodes is used.", "href": "../../../../source/kernels/PorousFlowMassTimeDerivative.html"}>>>
fluid_component<<<{"description": "The index corresponding to the component for this kernel"}>>> = 0
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pp
[]
[mass1]
type = PorousFlowMassTimeDerivative<<<{"description": "Derivative of fluid-component mass with respect to time. Mass lumping to the nodes is used.", "href": "../../../../source/kernels/PorousFlowMassTimeDerivative.html"}>>>
fluid_component<<<{"description": "The index corresponding to the component for this kernel"}>>> = 1
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = mass_frac_comp0
[]
[]
[UserObjects<<<{"href": "../../../../syntax/UserObjects/index.html"}>>>]
[dictator]
type = PorousFlowDictator<<<{"description": "Holds information on the PorousFlow variable names", "href": "../../../../source/userobjects/PorousFlowDictator.html"}>>>
porous_flow_vars<<<{"description": "List of primary variables that are used in the PorousFlow simulation. Jacobian entries involving derivatives wrt these variables will be computed. In single-phase models you will just have one (eg 'pressure'), in two-phase models you will have two (eg 'p_water p_gas', or 'p_water s_water'), etc."}>>> = 'pp mass_frac_comp0'
number_fluid_phases<<<{"description": "The number of fluid phases in the simulation"}>>> = 1
number_fluid_components<<<{"description": "The number of fluid components in the simulation"}>>> = 2
[]
[pc]
type = PorousFlowCapillaryPressureVG<<<{"description": "van Genuchten capillary pressure", "href": "../../../../source/userobjects/PorousFlowCapillaryPressureVG.html"}>>>
m<<<{"description": "van Genuchten exponent m. Must be between 0 and 1, and optimally should be set to >0.5"}>>> = 0.5
alpha<<<{"description": "van Genuchten parameter alpha. Must be positive"}>>> = 1
[]
[]
[FluidProperties<<<{"href": "../../../../syntax/FluidProperties/index.html"}>>>]
[simple_fluid]
type = SimpleFluidProperties<<<{"description": "Fluid properties for a simple fluid with a constant bulk density", "href": "../../../../source/fluidproperties/SimpleFluidProperties.html"}>>>
bulk_modulus<<<{"description": "Constant bulk modulus (Pa)"}>>> = 1
density0<<<{"description": "Density at zero pressure and zero temperature"}>>> = 1
thermal_expansion<<<{"description": "Constant coefficient of thermal expansion (1/K)"}>>> = 0
[]
[]
[Materials<<<{"href": "../../../../syntax/Materials/index.html"}>>>]
[temperature]
type = PorousFlowTemperature<<<{"description": "Material to provide temperature at the quadpoints or nodes and derivatives of it with respect to the PorousFlow variables", "href": "../../../../source/materials/PorousFlowTemperature.html"}>>>
[]
[ppss]
type = PorousFlow1PhaseP<<<{"description": "This Material is used for the partially saturated single-phase situation where porepressure is the primary variable", "href": "../../../../source/materials/PorousFlow1PhaseP.html"}>>>
porepressure<<<{"description": "Variable that represents the porepressure of the single phase"}>>> = pp
capillary_pressure<<<{"description": "Name of the UserObject defining the capillary pressure"}>>> = pc
[]
[massfrac]
type = PorousFlowMassFraction<<<{"description": "This Material forms a std::vector<std::vector ...> of mass-fractions out of the individual mass fractions", "href": "../../../../source/materials/PorousFlowMassFraction.html"}>>>
mass_fraction_vars<<<{"description": "List of variables that represent the mass fractions. Format is 'f_ph0^c0 f_ph0^c1 f_ph0^c2 ... f_ph0^c(N-2) f_ph1^c0 f_ph1^c1 fph1^c2 ... fph1^c(N-2) ... fphP^c0 f_phP^c1 fphP^c2 ... fphP^c(N-2)' where N=num_components and P=num_phases, and it is assumed that f_ph^c(N-1)=1-sum(f_ph^c,{c,0,N-2}) so that f_ph^c(N-1) need not be given. If no variables are provided then num_phases=1=num_components."}>>> = 'mass_frac_comp0'
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid<<<{"description": "This Material calculates fluid properties at the quadpoints or nodes for a single component fluid", "href": "../../../../source/materials/PorousFlowSingleComponentFluid.html"}>>>
fp<<<{"description": "The name of the user object for fluid properties"}>>> = simple_fluid
phase<<<{"description": "The phase number"}>>> = 0
[]
[porosity]
type = PorousFlowPorosityConst<<<{"description": "This Material calculates the porosity assuming it is constant", "href": "../../../../source/materials/PorousFlowPorosityConst.html"}>>>
porosity<<<{"description": "The porosity (assumed indepenent of porepressure, temperature, strain, etc, for this material). This should be a real number, or a constant monomial variable (not a linear lagrange or other kind of variable)."}>>> = 0.1
[]
[]
[Postprocessors<<<{"href": "../../../../syntax/Postprocessors/index.html"}>>>]
[total_mass_0]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
[]
[total_mass_1]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 1
[]
[]
[Preconditioning<<<{"href": "../../../../syntax/Preconditioning/index.html"}>>>]
[andy]
type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../../../../source/preconditioners/SingleMatrixPreconditioner.html"}>>>
full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
petsc_options_iname<<<{"description": "Names of PETSc name/value pairs"}>>> = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it'
petsc_options_value<<<{"description": "Values of PETSc name/value pairs (must correspond with \"petsc_options_iname\""}>>> = 'bcgs bjacobi 1 .999 10000'
[]
[]
[Executioner<<<{"href": "../../../../syntax/Executioner/index.html"}>>>]
type = Transient
solve_type = Newton
dt = 1
end_time = 1
[]
[Outputs<<<{"href": "../../../../syntax/Outputs/index.html"}>>>]
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'timestep_end'
file_base<<<{"description": "Common file base name to be utilized with all output objects"}>>> = mass02
csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
[]
(moose/modules/porous_flow/test/tests/mass_conservation/mass02.i)Table 3: Nodal masses in the single-phase, 2-component test
Density | Saturation | Nodal mass | Nodal mass | |||
---|---|---|---|---|---|---|
-1 | -1 | 0.367879441 | 0.707106781 | 1 | 0.008671 | 0 |
-0.333333333 | -0.333333333 | 0.716531311 | 0.948683298 | 0.111111 | 0.00251763 | 0.02014108 |
-0.333333333 | -0.333333333 | 0.716531311 | 0.948683298 | 0.111111 | 0.00251763 | 0.02014108 |
0.333333333 | 0.333333333 | 1.395612425 | 1 | 0.111111 | 0.00516893 | 0.04135148 |
0.333333333 | 0.333333333 | 1.395612425 | 1 | 0.111111 | 0.00516893 | 0.04135148 |
1 | 1 | 2.718281828 | 1 | 1 | 0.09060939 | 0 |
- | - | - | - | Total | 0.11465353 | 0.12298511 |
MOOSE produces the expected answer.
Two-phase, two-components
A 1D model with two elements from is created, with two phases (0 and 1), and two fluid components (0 and 1). The phase densities are calculated using a constant bulk modulus fluid with bulk modulus of 1 Pa. The density at zero pressure is 1 kg.m for phase 0, and 0.1 kg.m for phase 1. The porepressure of phase 0 is held fixed at 1 Pa, and a constant capillary pressure of 0 is specified so that the pressure of phase 1 is also 1 Pa. This results in phase densities of kg.m for phase 0, and kg.m for phase 1.
Saturation of phase 1 varies linearly as , while porosity is 0.1 throughout. The mass fraction of species 0 in fluid phase 0 is specified as 0.3, while the mass fraction of species 0 in phase 1 is 0.55.
# Checking that the mass postprocessor correctly calculates the mass
# of each component in each phase, as well as the total mass of each
# component in all phases.
# 2phase, 2component, constant porosity
[Mesh<<<{"href": "../../../../syntax/Mesh/index.html"}>>>]
type = GeneratedMesh
dim = 1
nx = 2
xmin = 0
xmax = 1
[]
[GlobalParams<<<{"href": "../../../../syntax/GlobalParams/index.html"}>>>]
PorousFlowDictator = dictator
[]
[Variables<<<{"href": "../../../../syntax/Variables/index.html"}>>>]
[pp]
[]
[sat]
[]
[]
[AuxVariables<<<{"href": "../../../../syntax/AuxVariables/index.html"}>>>]
[massfrac_ph0_sp0]
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 0.3
[]
[massfrac_ph1_sp0]
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 0.55
[]
[]
[ICs<<<{"href": "../../../../syntax/ICs/index.html"}>>>]
[pinit]
type = ConstantIC<<<{"description": "Sets a constant field value.", "href": "../../../../source/ics/ConstantIC.html"}>>>
value<<<{"description": "The value to be set in IC"}>>> = 1
variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = pp
[]
[satinit]
type = FunctionIC<<<{"description": "An initial condition that uses a normal function of x, y, z to produce values (and optionally gradients) for a field variable.", "href": "../../../../source/ics/FunctionIC.html"}>>>
function<<<{"description": "The initial condition function."}>>> = 1-x
variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = sat
[]
[]
[Kernels<<<{"href": "../../../../syntax/Kernels/index.html"}>>>]
[mass0]
type = PorousFlowMassTimeDerivative<<<{"description": "Derivative of fluid-component mass with respect to time. Mass lumping to the nodes is used.", "href": "../../../../source/kernels/PorousFlowMassTimeDerivative.html"}>>>
fluid_component<<<{"description": "The index corresponding to the component for this kernel"}>>> = 0
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pp
[]
[mass1]
type = PorousFlowMassTimeDerivative<<<{"description": "Derivative of fluid-component mass with respect to time. Mass lumping to the nodes is used.", "href": "../../../../source/kernels/PorousFlowMassTimeDerivative.html"}>>>
fluid_component<<<{"description": "The index corresponding to the component for this kernel"}>>> = 1
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = sat
[]
[]
[UserObjects<<<{"href": "../../../../syntax/UserObjects/index.html"}>>>]
[dictator]
type = PorousFlowDictator<<<{"description": "Holds information on the PorousFlow variable names", "href": "../../../../source/userobjects/PorousFlowDictator.html"}>>>
porous_flow_vars<<<{"description": "List of primary variables that are used in the PorousFlow simulation. Jacobian entries involving derivatives wrt these variables will be computed. In single-phase models you will just have one (eg 'pressure'), in two-phase models you will have two (eg 'p_water p_gas', or 'p_water s_water'), etc."}>>> = 'pp sat'
number_fluid_phases<<<{"description": "The number of fluid phases in the simulation"}>>> = 2
number_fluid_components<<<{"description": "The number of fluid components in the simulation"}>>> = 2
[]
[pc]
type = PorousFlowCapillaryPressureConst<<<{"description": "Constant capillary pressure", "href": "../../../../source/userobjects/PorousFlowCapillaryPressureConst.html"}>>>
pc<<<{"description": "Constant capillary pressure (Pa). Default is 0"}>>> = 0
[]
[]
[FluidProperties<<<{"href": "../../../../syntax/FluidProperties/index.html"}>>>]
[simple_fluid0]
type = SimpleFluidProperties<<<{"description": "Fluid properties for a simple fluid with a constant bulk density", "href": "../../../../source/fluidproperties/SimpleFluidProperties.html"}>>>
bulk_modulus<<<{"description": "Constant bulk modulus (Pa)"}>>> = 1
density0<<<{"description": "Density at zero pressure and zero temperature"}>>> = 1
thermal_expansion<<<{"description": "Constant coefficient of thermal expansion (1/K)"}>>> = 0
[]
[simple_fluid1]
type = SimpleFluidProperties<<<{"description": "Fluid properties for a simple fluid with a constant bulk density", "href": "../../../../source/fluidproperties/SimpleFluidProperties.html"}>>>
bulk_modulus<<<{"description": "Constant bulk modulus (Pa)"}>>> = 1
density0<<<{"description": "Density at zero pressure and zero temperature"}>>> = 0.1
thermal_expansion<<<{"description": "Constant coefficient of thermal expansion (1/K)"}>>> = 0
[]
[]
[Materials<<<{"href": "../../../../syntax/Materials/index.html"}>>>]
[temperature]
type = PorousFlowTemperature<<<{"description": "Material to provide temperature at the quadpoints or nodes and derivatives of it with respect to the PorousFlow variables", "href": "../../../../source/materials/PorousFlowTemperature.html"}>>>
[]
[ppss]
type = PorousFlow2PhasePS<<<{"description": "This Material calculates the 2 porepressures and the 2 saturations in a 2-phase situation, and derivatives of these with respect to the PorousFlowVariables.", "href": "../../../../source/materials/PorousFlow2PhasePS.html"}>>>
phase0_porepressure<<<{"description": "Variable that is the porepressure of phase 0 (the liquid phase)"}>>> = pp
phase1_saturation<<<{"description": "Variable that is the saturation of phase 1 (the gas phase)"}>>> = sat
capillary_pressure<<<{"description": "Name of the UserObject defining the capillary pressure"}>>> = pc
[]
[massfrac]
type = PorousFlowMassFraction<<<{"description": "This Material forms a std::vector<std::vector ...> of mass-fractions out of the individual mass fractions", "href": "../../../../source/materials/PorousFlowMassFraction.html"}>>>
mass_fraction_vars<<<{"description": "List of variables that represent the mass fractions. Format is 'f_ph0^c0 f_ph0^c1 f_ph0^c2 ... f_ph0^c(N-2) f_ph1^c0 f_ph1^c1 fph1^c2 ... fph1^c(N-2) ... fphP^c0 f_phP^c1 fphP^c2 ... fphP^c(N-2)' where N=num_components and P=num_phases, and it is assumed that f_ph^c(N-1)=1-sum(f_ph^c,{c,0,N-2}) so that f_ph^c(N-1) need not be given. If no variables are provided then num_phases=1=num_components."}>>> = 'massfrac_ph0_sp0 massfrac_ph1_sp0'
[]
[simple_fluid0]
type = PorousFlowSingleComponentFluid<<<{"description": "This Material calculates fluid properties at the quadpoints or nodes for a single component fluid", "href": "../../../../source/materials/PorousFlowSingleComponentFluid.html"}>>>
fp<<<{"description": "The name of the user object for fluid properties"}>>> = simple_fluid0
phase<<<{"description": "The phase number"}>>> = 0
[]
[simple_fluid1]
type = PorousFlowSingleComponentFluid<<<{"description": "This Material calculates fluid properties at the quadpoints or nodes for a single component fluid", "href": "../../../../source/materials/PorousFlowSingleComponentFluid.html"}>>>
fp<<<{"description": "The name of the user object for fluid properties"}>>> = simple_fluid1
phase<<<{"description": "The phase number"}>>> = 1
[]
[porosity]
type = PorousFlowPorosityConst<<<{"description": "This Material calculates the porosity assuming it is constant", "href": "../../../../source/materials/PorousFlowPorosityConst.html"}>>>
porosity<<<{"description": "The porosity (assumed indepenent of porepressure, temperature, strain, etc, for this material). This should be a real number, or a constant monomial variable (not a linear lagrange or other kind of variable)."}>>> = 0.1
[]
[]
[Postprocessors<<<{"href": "../../../../syntax/Postprocessors/index.html"}>>>]
[comp0_phase0_mass]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 0
phase<<<{"description": "The index of the fluid phase that this Postprocessor is restricted to. Multiple indices can be entered"}>>> = 0
[]
[comp0_phase1_mass]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 0
phase<<<{"description": "The index of the fluid phase that this Postprocessor is restricted to. Multiple indices can be entered"}>>> = 1
[]
[comp0_total_mass]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 0
[]
[comp0_total_mass2]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 0
phase<<<{"description": "The index of the fluid phase that this Postprocessor is restricted to. Multiple indices can be entered"}>>> = '0 1'
[]
[comp1_phase0_mass]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 1
phase<<<{"description": "The index of the fluid phase that this Postprocessor is restricted to. Multiple indices can be entered"}>>> = 0
[]
[comp1_phase1_mass]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 1
phase<<<{"description": "The index of the fluid phase that this Postprocessor is restricted to. Multiple indices can be entered"}>>> = 1
[]
[comp1_total_mass]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 1
[]
[comp1_total_mass2]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 1
phase<<<{"description": "The index of the fluid phase that this Postprocessor is restricted to. Multiple indices can be entered"}>>> = '0 1'
[]
[]
[Executioner<<<{"href": "../../../../syntax/Executioner/index.html"}>>>]
type = Transient
solve_type = Newton
nl_abs_tol = 1e-16
dt = 1
end_time = 1
[]
[Outputs<<<{"href": "../../../../syntax/Outputs/index.html"}>>>]
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'timestep_end'
file_base<<<{"description": "Common file base name to be utilized with all output objects"}>>> = mass05
csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
[]
(moose/modules/porous_flow/test/tests/mass_conservation/mass05.i)It is simple to calculate the total mass of each component in each phase using Eq. (1), the results of which are
Table 4: Masses in the 2-phase, 2-component test
Species | Phase | Total mass Eq. (1) | Total mass (MOOSE) |
---|---|---|---|
0 | 0 | 0.04077423 | 0.04077423 |
0 | 1 | 0.007475275 | 0.007475275 |
0 | all | 0.04824950 | 0.04824950 |
1 | 0 | 0.09513986 | 0.09513986 |
1 | 1 | 0.006116134 | 0.006116134 |
1 | all | 0.10125560 | 0.1012560 |
Constant fluid source
A fluid source of kg.m.s is introduced into a single element, and the PorousFlowFluidMass postprocessor is used to record the fluid mass as a function of time:
# checking that the mass postprocessor correctly calculates the mass
# 1phase, 1component, constant porosity, with a constant fluid source
[Mesh<<<{"href": "../../../../syntax/Mesh/index.html"}>>>]
type = GeneratedMesh
dim = 3
[]
[GlobalParams<<<{"href": "../../../../syntax/GlobalParams/index.html"}>>>]
PorousFlowDictator = dictator
[]
[Variables<<<{"href": "../../../../syntax/Variables/index.html"}>>>]
[pp]
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = -0.5
[]
[]
[Kernels<<<{"href": "../../../../syntax/Kernels/index.html"}>>>]
[mass0]
type = PorousFlowMassTimeDerivative<<<{"description": "Derivative of fluid-component mass with respect to time. Mass lumping to the nodes is used.", "href": "../../../../source/kernels/PorousFlowMassTimeDerivative.html"}>>>
fluid_component<<<{"description": "The index corresponding to the component for this kernel"}>>> = 0
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pp
[]
[source]
type = BodyForce<<<{"description": "Demonstrates the multiple ways that scalar values can be introduced into kernels, e.g. (controllable) constants, functions, and postprocessors. Implements the weak form $(\\psi_i, -f)$.", "href": "../../../../source/kernels/BodyForce.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pp
value<<<{"description": "Coefficient to multiply by the body force term"}>>> = 0.1 # kg/m^3/s
[]
[]
[UserObjects<<<{"href": "../../../../syntax/UserObjects/index.html"}>>>]
[dictator]
type = PorousFlowDictator<<<{"description": "Holds information on the PorousFlow variable names", "href": "../../../../source/userobjects/PorousFlowDictator.html"}>>>
porous_flow_vars<<<{"description": "List of primary variables that are used in the PorousFlow simulation. Jacobian entries involving derivatives wrt these variables will be computed. In single-phase models you will just have one (eg 'pressure'), in two-phase models you will have two (eg 'p_water p_gas', or 'p_water s_water'), etc."}>>> = 'pp'
number_fluid_phases<<<{"description": "The number of fluid phases in the simulation"}>>> = 1
number_fluid_components<<<{"description": "The number of fluid components in the simulation"}>>> = 1
[]
[pc]
type = PorousFlowCapillaryPressureVG<<<{"description": "van Genuchten capillary pressure", "href": "../../../../source/userobjects/PorousFlowCapillaryPressureVG.html"}>>>
m<<<{"description": "van Genuchten exponent m. Must be between 0 and 1, and optimally should be set to >0.5"}>>> = 0.5
alpha<<<{"description": "van Genuchten parameter alpha. Must be positive"}>>> = 1
[]
[]
[FluidProperties<<<{"href": "../../../../syntax/FluidProperties/index.html"}>>>]
[simple_fluid]
type = SimpleFluidProperties<<<{"description": "Fluid properties for a simple fluid with a constant bulk density", "href": "../../../../source/fluidproperties/SimpleFluidProperties.html"}>>>
bulk_modulus<<<{"description": "Constant bulk modulus (Pa)"}>>> = 1
density0<<<{"description": "Density at zero pressure and zero temperature"}>>> = 1
thermal_expansion<<<{"description": "Constant coefficient of thermal expansion (1/K)"}>>> = 0
[]
[]
[Materials<<<{"href": "../../../../syntax/Materials/index.html"}>>>]
[temperature]
type = PorousFlowTemperature<<<{"description": "Material to provide temperature at the quadpoints or nodes and derivatives of it with respect to the PorousFlow variables", "href": "../../../../source/materials/PorousFlowTemperature.html"}>>>
[]
[ppss]
type = PorousFlow1PhaseP<<<{"description": "This Material is used for the partially saturated single-phase situation where porepressure is the primary variable", "href": "../../../../source/materials/PorousFlow1PhaseP.html"}>>>
porepressure<<<{"description": "Variable that represents the porepressure of the single phase"}>>> = pp
capillary_pressure<<<{"description": "Name of the UserObject defining the capillary pressure"}>>> = pc
[]
[massfrac]
type = PorousFlowMassFraction<<<{"description": "This Material forms a std::vector<std::vector ...> of mass-fractions out of the individual mass fractions", "href": "../../../../source/materials/PorousFlowMassFraction.html"}>>>
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid<<<{"description": "This Material calculates fluid properties at the quadpoints or nodes for a single component fluid", "href": "../../../../source/materials/PorousFlowSingleComponentFluid.html"}>>>
fp<<<{"description": "The name of the user object for fluid properties"}>>> = simple_fluid
phase<<<{"description": "The phase number"}>>> = 0
[]
[porosity]
type = PorousFlowPorosityConst<<<{"description": "This Material calculates the porosity assuming it is constant", "href": "../../../../source/materials/PorousFlowPorosityConst.html"}>>>
porosity<<<{"description": "The porosity (assumed indepenent of porepressure, temperature, strain, etc, for this material). This should be a real number, or a constant monomial variable (not a linear lagrange or other kind of variable)."}>>> = 0.1
[]
[]
[Postprocessors<<<{"href": "../../../../syntax/Postprocessors/index.html"}>>>]
[porepressure]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../../source/postprocessors/PointValue.html"}>>>
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0 0 0'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = pp
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
[]
[total_mass]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
[]
[]
[Preconditioning<<<{"href": "../../../../syntax/Preconditioning/index.html"}>>>]
[andy]
type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../../../../source/preconditioners/SingleMatrixPreconditioner.html"}>>>
full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
petsc_options_iname<<<{"description": "Names of PETSc name/value pairs"}>>> = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it'
petsc_options_value<<<{"description": "Values of PETSc name/value pairs (must correspond with \"petsc_options_iname\""}>>> = 'gmres bjacobi 1E-12 1E-20 10000'
[]
[]
[Executioner<<<{"href": "../../../../syntax/Executioner/index.html"}>>>]
type = Transient
solve_type = Newton
dt = 1
end_time = 10
[]
[Outputs<<<{"href": "../../../../syntax/Outputs/index.html"}>>>]
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
file_base<<<{"description": "Common file base name to be utilized with all output objects"}>>> = mass03
csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
[]
(moose/modules/porous_flow/test/tests/mass_conservation/mass03.i)Mass conservation in a deforming material
A single unit element, with roller BCs on its sides and bottom, is compressed at a uniform rate: Here is the vertical displacement and is time. There is no fluid flow and the material's boundaries are impermeable. Fluid mass conservation is checked.
# The sample is a single unit element, with roller BCs on the sides
# and bottom. A constant displacement is applied to the top: disp_z = -0.01*t.
# There is no fluid flow.
# Fluid mass conservation is checked.
#
# Under these conditions
# porepressure = porepressure(t=0) - (Fluid bulk modulus)*log(1 - 0.01*t)
# stress_xx = (bulk - 2*shear/3)*disp_z/L (remember this is effective stress)
# stress_zz = (bulk + 4*shear/3)*disp_z/L (remember this is effective stress)
# where L is the height of the sample (L=1 in this test)
#
# Parameters:
# Bulk modulus = 2
# Shear modulus = 1.5
# fluid bulk modulus = 0.5
# initial porepressure = 0.1
#
# Desired output:
# zdisp = -0.01*t
# p0 = 0.1 - 0.5*log(1-0.01*t)
# stress_xx = stress_yy = -0.01*t
# stress_zz = -0.04*t
#
# Regarding the "log" - it comes from preserving fluid mass
[Mesh<<<{"href": "../../../../syntax/Mesh/index.html"}>>>]
type = GeneratedMesh
dim = 3
nx = 1
ny = 1
nz = 1
xmin = -0.5
xmax = 0.5
ymin = -0.5
ymax = 0.5
zmin = -0.5
zmax = 0.5
[]
[GlobalParams<<<{"href": "../../../../syntax/GlobalParams/index.html"}>>>]
displacements = 'disp_x disp_y disp_z'
PorousFlowDictator = dictator
block = 0
[]
[Variables<<<{"href": "../../../../syntax/Variables/index.html"}>>>]
[disp_x]
[]
[disp_y]
[]
[disp_z]
[]
[porepressure]
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 0.1
[]
[]
[BCs<<<{"href": "../../../../syntax/BCs/index.html"}>>>]
[confinex]
type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../../source/bcs/DirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
value<<<{"description": "Value of the BC"}>>> = 0
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left right'
[]
[confiney]
type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../../source/bcs/DirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
value<<<{"description": "Value of the BC"}>>> = 0
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom top'
[]
[basefixed]
type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../../source/bcs/DirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
value<<<{"description": "Value of the BC"}>>> = 0
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = back
[]
[top_velocity]
type = FunctionDirichletBC<<<{"description": "Imposes the essential boundary condition $u=g(t,\\vec{x})$, where $g$ is a (possibly) time and space-dependent MOOSE Function.", "href": "../../../../source/bcs/FunctionDirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
function<<<{"description": "The forcing function."}>>> = -0.01*t
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = front
[]
[]
[Kernels<<<{"href": "../../../../syntax/Kernels/index.html"}>>>]
[grad_stress_x]
type = StressDivergenceTensors<<<{"description": "Stress divergence kernel for the Cartesian coordinate system", "href": "../../../../source/kernels/StressDivergenceTensors.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
component<<<{"description": "An integer corresponding to the direction the variable this kernel acts in. (0 for x, 1 for y, 2 for z)"}>>> = 0
[]
[grad_stress_y]
type = StressDivergenceTensors<<<{"description": "Stress divergence kernel for the Cartesian coordinate system", "href": "../../../../source/kernels/StressDivergenceTensors.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
component<<<{"description": "An integer corresponding to the direction the variable this kernel acts in. (0 for x, 1 for y, 2 for z)"}>>> = 1
[]
[grad_stress_z]
type = StressDivergenceTensors<<<{"description": "Stress divergence kernel for the Cartesian coordinate system", "href": "../../../../source/kernels/StressDivergenceTensors.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
component<<<{"description": "An integer corresponding to the direction the variable this kernel acts in. (0 for x, 1 for y, 2 for z)"}>>> = 2
[]
[poro_x]
type = PorousFlowEffectiveStressCoupling<<<{"description": "Implements the weak form of the expression biot_coefficient * grad(effective fluid pressure)", "href": "../../../../source/kernels/PorousFlowEffectiveStressCoupling.html"}>>>
biot_coefficient<<<{"description": "Biot coefficient"}>>> = 0.3
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
component<<<{"description": "The component (0 for x, 1 for y and 2 for z) of grad(P)"}>>> = 0
[]
[poro_y]
type = PorousFlowEffectiveStressCoupling<<<{"description": "Implements the weak form of the expression biot_coefficient * grad(effective fluid pressure)", "href": "../../../../source/kernels/PorousFlowEffectiveStressCoupling.html"}>>>
biot_coefficient<<<{"description": "Biot coefficient"}>>> = 0.3
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
component<<<{"description": "The component (0 for x, 1 for y and 2 for z) of grad(P)"}>>> = 1
[]
[poro_z]
type = PorousFlowEffectiveStressCoupling<<<{"description": "Implements the weak form of the expression biot_coefficient * grad(effective fluid pressure)", "href": "../../../../source/kernels/PorousFlowEffectiveStressCoupling.html"}>>>
biot_coefficient<<<{"description": "Biot coefficient"}>>> = 0.3
component<<<{"description": "The component (0 for x, 1 for y and 2 for z) of grad(P)"}>>> = 2
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
[]
[poro_vol_exp]
type = PorousFlowMassVolumetricExpansion<<<{"description": "Component_mass*rate_of_solid_volumetric_expansion. This Kernel lumps the component mass to the nodes.", "href": "../../../../source/kernels/PorousFlowMassVolumetricExpansion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = porepressure
fluid_component<<<{"description": "The index corresponding to the component for this kernel"}>>> = 0
[]
[mass0]
type = PorousFlowMassTimeDerivative<<<{"description": "Derivative of fluid-component mass with respect to time. Mass lumping to the nodes is used.", "href": "../../../../source/kernels/PorousFlowMassTimeDerivative.html"}>>>
fluid_component<<<{"description": "The index corresponding to the component for this kernel"}>>> = 0
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = porepressure
[]
[]
[AuxVariables<<<{"href": "../../../../syntax/AuxVariables/index.html"}>>>]
[stress_xx]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[]
[stress_xy]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[]
[stress_xz]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[]
[stress_yy]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[]
[stress_yz]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[]
[stress_zz]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[]
[]
[AuxKernels<<<{"href": "../../../../syntax/AuxKernels/index.html"}>>>]
[stress_xx]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../../../../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_xx
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 0
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 0
[]
[stress_xy]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../../../../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_xy
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 0
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 1
[]
[stress_xz]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../../../../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_xz
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 0
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 2
[]
[stress_yy]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../../../../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_yy
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 1
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 1
[]
[stress_yz]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../../../../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_yz
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 1
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 2
[]
[stress_zz]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../../../../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_zz
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 2
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 2
[]
[]
[FluidProperties<<<{"href": "../../../../syntax/FluidProperties/index.html"}>>>]
[simple_fluid]
type = SimpleFluidProperties<<<{"description": "Fluid properties for a simple fluid with a constant bulk density", "href": "../../../../source/fluidproperties/SimpleFluidProperties.html"}>>>
bulk_modulus<<<{"description": "Constant bulk modulus (Pa)"}>>> = 0.5
density0<<<{"description": "Density at zero pressure and zero temperature"}>>> = 1
thermal_expansion<<<{"description": "Constant coefficient of thermal expansion (1/K)"}>>> = 0
viscosity<<<{"description": "Constant dynamic viscosity (Pa.s)"}>>> = 1
[]
[]
[Materials<<<{"href": "../../../../syntax/Materials/index.html"}>>>]
[temperature]
type = PorousFlowTemperature<<<{"description": "Material to provide temperature at the quadpoints or nodes and derivatives of it with respect to the PorousFlow variables", "href": "../../../../source/materials/PorousFlowTemperature.html"}>>>
[]
[elasticity_tensor]
type = ComputeElasticityTensor<<<{"description": "Compute an elasticity tensor.", "href": "../../../../source/materials/ComputeElasticityTensor.html"}>>>
C_ijkl<<<{"description": "Stiffness tensor for material"}>>> = '1 1.5'
# bulk modulus is lambda + 2*mu/3 = 1 + 2*1.5/3 = 2
fill_method<<<{"description": "The fill method"}>>> = symmetric_isotropic
[]
[strain]
type = ComputeSmallStrain<<<{"description": "Compute a small strain.", "href": "../../../../source/materials/ComputeSmallStrain.html"}>>>
[]
[stress]
type = ComputeLinearElasticStress<<<{"description": "Compute stress using elasticity for small strains", "href": "../../../../source/materials/ComputeLinearElasticStress.html"}>>>
[]
[vol_strain]
type = PorousFlowVolumetricStrain<<<{"description": "Compute volumetric strain and the volumetric_strain rate, for use in PorousFlow.", "href": "../../../../source/materials/PorousFlowVolumetricStrain.html"}>>>
[]
[eff_fluid_pressure]
type = PorousFlowEffectiveFluidPressure<<<{"description": "This Material calculates an effective fluid pressure: effective_stress = total_stress + biot_coeff*effective_fluid_pressure. The effective_fluid_pressure = sum_{phases}(S_phase * P_phase)", "href": "../../../../source/materials/PorousFlowEffectiveFluidPressure.html"}>>>
[]
[ppss]
type = PorousFlow1PhaseP<<<{"description": "This Material is used for the partially saturated single-phase situation where porepressure is the primary variable", "href": "../../../../source/materials/PorousFlow1PhaseP.html"}>>>
porepressure<<<{"description": "Variable that represents the porepressure of the single phase"}>>> = porepressure
capillary_pressure<<<{"description": "Name of the UserObject defining the capillary pressure"}>>> = pc
[]
[massfrac]
type = PorousFlowMassFraction<<<{"description": "This Material forms a std::vector<std::vector ...> of mass-fractions out of the individual mass fractions", "href": "../../../../source/materials/PorousFlowMassFraction.html"}>>>
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid<<<{"description": "This Material calculates fluid properties at the quadpoints or nodes for a single component fluid", "href": "../../../../source/materials/PorousFlowSingleComponentFluid.html"}>>>
fp<<<{"description": "The name of the user object for fluid properties"}>>> = simple_fluid
phase<<<{"description": "The phase number"}>>> = 0
[]
[porosity]
type = PorousFlowPorosityConst<<<{"description": "This Material calculates the porosity assuming it is constant", "href": "../../../../source/materials/PorousFlowPorosityConst.html"}>>>
porosity<<<{"description": "The porosity (assumed indepenent of porepressure, temperature, strain, etc, for this material). This should be a real number, or a constant monomial variable (not a linear lagrange or other kind of variable)."}>>> = 0.1
[]
[permeability]
type = PorousFlowPermeabilityConst<<<{"description": "This Material calculates the permeability tensor assuming it is constant", "href": "../../../../source/materials/PorousFlowPermeabilityConst.html"}>>>
permeability<<<{"description": "The permeability tensor (usually in m^2), which is assumed constant for this material"}>>> = '0.5 0 0 0 0.5 0 0 0 0.5'
[]
[]
[UserObjects<<<{"href": "../../../../syntax/UserObjects/index.html"}>>>]
[dictator]
type = PorousFlowDictator<<<{"description": "Holds information on the PorousFlow variable names", "href": "../../../../source/userobjects/PorousFlowDictator.html"}>>>
porous_flow_vars<<<{"description": "List of primary variables that are used in the PorousFlow simulation. Jacobian entries involving derivatives wrt these variables will be computed. In single-phase models you will just have one (eg 'pressure'), in two-phase models you will have two (eg 'p_water p_gas', or 'p_water s_water'), etc."}>>> = 'porepressure disp_x disp_y disp_z'
number_fluid_phases<<<{"description": "The number of fluid phases in the simulation"}>>> = 1
number_fluid_components<<<{"description": "The number of fluid components in the simulation"}>>> = 1
[]
[pc]
type = PorousFlowCapillaryPressureVG<<<{"description": "van Genuchten capillary pressure", "href": "../../../../source/userobjects/PorousFlowCapillaryPressureVG.html"}>>>
m<<<{"description": "van Genuchten exponent m. Must be between 0 and 1, and optimally should be set to >0.5"}>>> = 0.5
alpha<<<{"description": "van Genuchten parameter alpha. Must be positive"}>>> = 1
[]
[]
[Postprocessors<<<{"href": "../../../../syntax/Postprocessors/index.html"}>>>]
[p0]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../../source/postprocessors/PointValue.html"}>>>
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = 'console csv'
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0 0 0'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = porepressure
[]
[zdisp]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../../source/postprocessors/PointValue.html"}>>>
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0 0 0.5'
use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = false
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = disp_z
[]
[stress_xx]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../../source/postprocessors/PointValue.html"}>>>
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0 0 0'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = stress_xx
[]
[stress_yy]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../../source/postprocessors/PointValue.html"}>>>
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0 0 0'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = stress_yy
[]
[stress_zz]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../../source/postprocessors/PointValue.html"}>>>
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0 0 0'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = stress_zz
[]
[fluid_mass]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 0
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = 'console csv'
[]
[]
[Preconditioning<<<{"href": "../../../../syntax/Preconditioning/index.html"}>>>]
[andy]
type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../../../../source/preconditioners/SingleMatrixPreconditioner.html"}>>>
full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
petsc_options_iname<<<{"description": "Names of PETSc name/value pairs"}>>> = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it'
petsc_options_value<<<{"description": "Values of PETSc name/value pairs (must correspond with \"petsc_options_iname\""}>>> = 'bcgs bjacobi 1E-14 1E-8 10000'
[]
[]
[Executioner<<<{"href": "../../../../syntax/Executioner/index.html"}>>>]
type = Transient
solve_type = Newton
start_time = 0
end_time = 10
dt = 2
[]
[Outputs<<<{"href": "../../../../syntax/Outputs/index.html"}>>>]
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
file_base<<<{"description": "Common file base name to be utilized with all output objects"}>>> = mass04
[csv]
type = CSV<<<{"description": "Output for postprocessors, vector postprocessors, and scalar variables using comma seperated values (CSV).", "href": "../../../../source/outputs/CSV.html"}>>>
[]
[]
(moose/modules/porous_flow/test/tests/mass_conservation/mass04.i)Note the use_displaced_mesh = true
option in the PorousFlowFluidMass postprocessor: it is necessary to correctly compute the mass.
Under these conditions Here is the fluid porepressure, which is at ; is the fluid bulk modulus; is the effective stress; is the drained bulk modulus of the material; is the shear modulus of the material; and is the height of the sample.
PorousFlow produces these results exactly, and, importantly, conserves fluid mass.
Similar tests are run in "RZ" coordinates.
Mass conservation in a deformable material with a source
A single unit element, with roller BCs on its sides and bottom, is injected with fluid at rate kg.s. Its top is free to move.
# The sample is a single unit element, with roller BCs on the sides and bottom.
# The top is free to move and fluid is injected at a constant rate of 1kg/s
# There is no fluid flow.
# Fluid mass conservation is checked.
# Under these conditions the fluid mass should increase at 1kg/s
# The porepressure should increase: rho0 * exp(P/bulk) = rho * exp(P0/bulk) + 1*t
# The stress_zz should be exactly biot * P since total stress is zero
[Mesh<<<{"href": "../../../../syntax/Mesh/index.html"}>>>]
type = GeneratedMesh
dim = 3
nx = 1
ny = 1
nz = 1
xmin = -0.5
xmax = 0.5
ymin = -0.5
ymax = 0.5
zmin = -0.5
zmax = 0.5
[]
[GlobalParams<<<{"href": "../../../../syntax/GlobalParams/index.html"}>>>]
displacements = 'disp_x disp_y disp_z'
PorousFlowDictator = dictator
block = 0
[]
[Variables<<<{"href": "../../../../syntax/Variables/index.html"}>>>]
[disp_x]
[]
[disp_y]
[]
[disp_z]
[]
[porepressure]
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 0.1
[]
[]
[BCs<<<{"href": "../../../../syntax/BCs/index.html"}>>>]
[confinex]
type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../../source/bcs/DirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
value<<<{"description": "Value of the BC"}>>> = 0
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left right'
[]
[confiney]
type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../../source/bcs/DirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
value<<<{"description": "Value of the BC"}>>> = 0
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom top'
[]
[basefixed]
type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../../source/bcs/DirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
value<<<{"description": "Value of the BC"}>>> = 0
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = back
[]
[]
[Kernels<<<{"href": "../../../../syntax/Kernels/index.html"}>>>]
[grad_stress_x]
type = StressDivergenceTensors<<<{"description": "Stress divergence kernel for the Cartesian coordinate system", "href": "../../../../source/kernels/StressDivergenceTensors.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
component<<<{"description": "An integer corresponding to the direction the variable this kernel acts in. (0 for x, 1 for y, 2 for z)"}>>> = 0
[]
[grad_stress_y]
type = StressDivergenceTensors<<<{"description": "Stress divergence kernel for the Cartesian coordinate system", "href": "../../../../source/kernels/StressDivergenceTensors.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
component<<<{"description": "An integer corresponding to the direction the variable this kernel acts in. (0 for x, 1 for y, 2 for z)"}>>> = 1
[]
[grad_stress_z]
type = StressDivergenceTensors<<<{"description": "Stress divergence kernel for the Cartesian coordinate system", "href": "../../../../source/kernels/StressDivergenceTensors.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
component<<<{"description": "An integer corresponding to the direction the variable this kernel acts in. (0 for x, 1 for y, 2 for z)"}>>> = 2
[]
[poro_x]
type = PorousFlowEffectiveStressCoupling<<<{"description": "Implements the weak form of the expression biot_coefficient * grad(effective fluid pressure)", "href": "../../../../source/kernels/PorousFlowEffectiveStressCoupling.html"}>>>
biot_coefficient<<<{"description": "Biot coefficient"}>>> = 0.3
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
component<<<{"description": "The component (0 for x, 1 for y and 2 for z) of grad(P)"}>>> = 0
[]
[poro_y]
type = PorousFlowEffectiveStressCoupling<<<{"description": "Implements the weak form of the expression biot_coefficient * grad(effective fluid pressure)", "href": "../../../../source/kernels/PorousFlowEffectiveStressCoupling.html"}>>>
biot_coefficient<<<{"description": "Biot coefficient"}>>> = 0.3
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
component<<<{"description": "The component (0 for x, 1 for y and 2 for z) of grad(P)"}>>> = 1
[]
[poro_z]
type = PorousFlowEffectiveStressCoupling<<<{"description": "Implements the weak form of the expression biot_coefficient * grad(effective fluid pressure)", "href": "../../../../source/kernels/PorousFlowEffectiveStressCoupling.html"}>>>
biot_coefficient<<<{"description": "Biot coefficient"}>>> = 0.3
component<<<{"description": "The component (0 for x, 1 for y and 2 for z) of grad(P)"}>>> = 2
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
[]
[poro_vol_exp]
type = PorousFlowMassVolumetricExpansion<<<{"description": "Component_mass*rate_of_solid_volumetric_expansion. This Kernel lumps the component mass to the nodes.", "href": "../../../../source/kernels/PorousFlowMassVolumetricExpansion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = porepressure
fluid_component<<<{"description": "The index corresponding to the component for this kernel"}>>> = 0
[]
[mass0]
type = PorousFlowMassTimeDerivative<<<{"description": "Derivative of fluid-component mass with respect to time. Mass lumping to the nodes is used.", "href": "../../../../source/kernels/PorousFlowMassTimeDerivative.html"}>>>
fluid_component<<<{"description": "The index corresponding to the component for this kernel"}>>> = 0
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = porepressure
[]
[]
[DiracKernels<<<{"href": "../../../../syntax/DiracKernels/index.html"}>>>]
[inject]
type = PorousFlowPointSourceFromPostprocessor<<<{"description": "Point source (or sink) that adds (or removes) fluid at a mass flux rate specified by a postprocessor.", "href": "../../../../source/dirackernels/PorousFlowPointSourceFromPostprocessor.html"}>>>
point<<<{"description": "The x,y,z coordinates of the point source (or sink)"}>>> = '0 0 0'
mass_flux<<<{"description": "The postprocessor name holding the mass flux at this point in kg/s (positive is flux in, negative is flux out)"}>>> = 1.0
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = porepressure
[]
[]
[AuxVariables<<<{"href": "../../../../syntax/AuxVariables/index.html"}>>>]
[stress_xx]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[]
[stress_xy]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[]
[stress_xz]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[]
[stress_yy]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[]
[stress_yz]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[]
[stress_zz]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[]
[]
[AuxKernels<<<{"href": "../../../../syntax/AuxKernels/index.html"}>>>]
[stress_xx]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../../../../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_xx
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 0
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 0
[]
[stress_xy]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../../../../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_xy
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 0
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 1
[]
[stress_xz]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../../../../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_xz
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 0
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 2
[]
[stress_yy]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../../../../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_yy
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 1
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 1
[]
[stress_yz]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../../../../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_yz
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 1
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 2
[]
[stress_zz]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../../../../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_zz
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 2
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 2
[]
[]
[FluidProperties<<<{"href": "../../../../syntax/FluidProperties/index.html"}>>>]
[simple_fluid]
type = SimpleFluidProperties<<<{"description": "Fluid properties for a simple fluid with a constant bulk density", "href": "../../../../source/fluidproperties/SimpleFluidProperties.html"}>>>
bulk_modulus<<<{"description": "Constant bulk modulus (Pa)"}>>> = 0.5
density0<<<{"description": "Density at zero pressure and zero temperature"}>>> = 1
thermal_expansion<<<{"description": "Constant coefficient of thermal expansion (1/K)"}>>> = 0
viscosity<<<{"description": "Constant dynamic viscosity (Pa.s)"}>>> = 1
[]
[]
[Materials<<<{"href": "../../../../syntax/Materials/index.html"}>>>]
[temperature]
type = PorousFlowTemperature<<<{"description": "Material to provide temperature at the quadpoints or nodes and derivatives of it with respect to the PorousFlow variables", "href": "../../../../source/materials/PorousFlowTemperature.html"}>>>
[]
[elasticity_tensor]
type = ComputeElasticityTensor<<<{"description": "Compute an elasticity tensor.", "href": "../../../../source/materials/ComputeElasticityTensor.html"}>>>
C_ijkl<<<{"description": "Stiffness tensor for material"}>>> = '1 1.5'
# bulk modulus is lambda + 2*mu/3 = 1 + 2*1.5/3 = 2
fill_method<<<{"description": "The fill method"}>>> = symmetric_isotropic
[]
[strain]
type = ComputeSmallStrain<<<{"description": "Compute a small strain.", "href": "../../../../source/materials/ComputeSmallStrain.html"}>>>
[]
[stress]
type = ComputeLinearElasticStress<<<{"description": "Compute stress using elasticity for small strains", "href": "../../../../source/materials/ComputeLinearElasticStress.html"}>>>
[]
[vol_strain]
type = PorousFlowVolumetricStrain<<<{"description": "Compute volumetric strain and the volumetric_strain rate, for use in PorousFlow.", "href": "../../../../source/materials/PorousFlowVolumetricStrain.html"}>>>
[]
[eff_fluid_pressure]
type = PorousFlowEffectiveFluidPressure<<<{"description": "This Material calculates an effective fluid pressure: effective_stress = total_stress + biot_coeff*effective_fluid_pressure. The effective_fluid_pressure = sum_{phases}(S_phase * P_phase)", "href": "../../../../source/materials/PorousFlowEffectiveFluidPressure.html"}>>>
[]
[ppss]
type = PorousFlow1PhaseP<<<{"description": "This Material is used for the partially saturated single-phase situation where porepressure is the primary variable", "href": "../../../../source/materials/PorousFlow1PhaseP.html"}>>>
porepressure<<<{"description": "Variable that represents the porepressure of the single phase"}>>> = porepressure
capillary_pressure<<<{"description": "Name of the UserObject defining the capillary pressure"}>>> = pc
[]
[massfrac]
type = PorousFlowMassFraction<<<{"description": "This Material forms a std::vector<std::vector ...> of mass-fractions out of the individual mass fractions", "href": "../../../../source/materials/PorousFlowMassFraction.html"}>>>
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid<<<{"description": "This Material calculates fluid properties at the quadpoints or nodes for a single component fluid", "href": "../../../../source/materials/PorousFlowSingleComponentFluid.html"}>>>
fp<<<{"description": "The name of the user object for fluid properties"}>>> = simple_fluid
phase<<<{"description": "The phase number"}>>> = 0
[]
[porosity]
type = PorousFlowPorosityConst<<<{"description": "This Material calculates the porosity assuming it is constant", "href": "../../../../source/materials/PorousFlowPorosityConst.html"}>>>
porosity<<<{"description": "The porosity (assumed indepenent of porepressure, temperature, strain, etc, for this material). This should be a real number, or a constant monomial variable (not a linear lagrange or other kind of variable)."}>>> = 0.1
[]
[permeability]
type = PorousFlowPermeabilityConst<<<{"description": "This Material calculates the permeability tensor assuming it is constant", "href": "../../../../source/materials/PorousFlowPermeabilityConst.html"}>>>
permeability<<<{"description": "The permeability tensor (usually in m^2), which is assumed constant for this material"}>>> = '0.5 0 0 0 0.5 0 0 0 0.5'
[]
[]
[UserObjects<<<{"href": "../../../../syntax/UserObjects/index.html"}>>>]
[dictator]
type = PorousFlowDictator<<<{"description": "Holds information on the PorousFlow variable names", "href": "../../../../source/userobjects/PorousFlowDictator.html"}>>>
porous_flow_vars<<<{"description": "List of primary variables that are used in the PorousFlow simulation. Jacobian entries involving derivatives wrt these variables will be computed. In single-phase models you will just have one (eg 'pressure'), in two-phase models you will have two (eg 'p_water p_gas', or 'p_water s_water'), etc."}>>> = 'porepressure disp_x disp_y disp_z'
number_fluid_phases<<<{"description": "The number of fluid phases in the simulation"}>>> = 1
number_fluid_components<<<{"description": "The number of fluid components in the simulation"}>>> = 1
[]
[pc]
type = PorousFlowCapillaryPressureVG<<<{"description": "van Genuchten capillary pressure", "href": "../../../../source/userobjects/PorousFlowCapillaryPressureVG.html"}>>>
m<<<{"description": "van Genuchten exponent m. Must be between 0 and 1, and optimally should be set to >0.5"}>>> = 0.5
alpha<<<{"description": "van Genuchten parameter alpha. Must be positive"}>>> = 1
[]
[]
[Postprocessors<<<{"href": "../../../../syntax/Postprocessors/index.html"}>>>]
[p0]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../../source/postprocessors/PointValue.html"}>>>
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = 'console csv'
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0 0 0'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = porepressure
[]
[zdisp]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../../source/postprocessors/PointValue.html"}>>>
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0 0 0.5'
use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = false
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = disp_z
[]
[stress_xx]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../../source/postprocessors/PointValue.html"}>>>
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0 0 0'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = stress_xx
[]
[stress_yy]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../../source/postprocessors/PointValue.html"}>>>
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0 0 0'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = stress_yy
[]
[stress_zz]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../../source/postprocessors/PointValue.html"}>>>
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0 0 0'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = stress_zz
[]
[fluid_mass]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 0
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = 'console csv'
[]
[]
[Preconditioning<<<{"href": "../../../../syntax/Preconditioning/index.html"}>>>]
[andy]
type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../../../../source/preconditioners/SingleMatrixPreconditioner.html"}>>>
full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
petsc_options_iname<<<{"description": "Names of PETSc name/value pairs"}>>> = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it'
petsc_options_value<<<{"description": "Values of PETSc name/value pairs (must correspond with \"petsc_options_iname\""}>>> = 'bcgs bjacobi 1E-14 1E-8 10000'
[]
[]
[Executioner<<<{"href": "../../../../syntax/Executioner/index.html"}>>>]
type = Transient
solve_type = Newton
start_time = 0
end_time = 10
dt = 2
[]
[Outputs<<<{"href": "../../../../syntax/Outputs/index.html"}>>>]
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
[csv]
type = CSV<<<{"description": "Output for postprocessors, vector postprocessors, and scalar variables using comma seperated values (CSV).", "href": "../../../../source/outputs/CSV.html"}>>>
[]
[]
(moose/modules/porous_flow/test/tests/mass_conservation/mass11.i)Under these conditions the fluid mass should increase at rate , and the porepressure should increase accordingly. The total stress should be zero since the top is free to move, so the effective stress should be . MOOSE produces these results exactly.
Similar tests are run in "RZ" coordinates.
Mass computation with a saturation threshold in multi-component, multi-phase fluids
The PorousFlowFluidMass postprocessor may be used to compute the total mass of each component in each phase, as well as the total mass of each component in all phases. Furthermore, a saturation threshold may be set to only count the fluid above the threshold.
# Checking that the mass postprocessor correctly calculates the mass
# of each component in each phase, as well as the total mass of each
# component in all phases. Also tests that optional saturation threshold
# gives the correct mass
# 2phase, 2component, constant porosity
# saturation_threshold set to 0.6 for phase 1
[Mesh<<<{"href": "../../../../syntax/Mesh/index.html"}>>>]
type = GeneratedMesh
dim = 1
nx = 10
xmin = 0
xmax = 1
[]
[GlobalParams<<<{"href": "../../../../syntax/GlobalParams/index.html"}>>>]
PorousFlowDictator = dictator
[]
[Variables<<<{"href": "../../../../syntax/Variables/index.html"}>>>]
[pp]
[]
[sat]
[]
[]
[AuxVariables<<<{"href": "../../../../syntax/AuxVariables/index.html"}>>>]
[massfrac_ph0_sp0]
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1
[]
[massfrac_ph1_sp0]
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 0
[]
[]
[ICs<<<{"href": "../../../../syntax/ICs/index.html"}>>>]
[pinit]
type = ConstantIC<<<{"description": "Sets a constant field value.", "href": "../../../../source/ics/ConstantIC.html"}>>>
value<<<{"description": "The value to be set in IC"}>>> = 1
variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = pp
[]
[satinit]
type = FunctionIC<<<{"description": "An initial condition that uses a normal function of x, y, z to produce values (and optionally gradients) for a field variable.", "href": "../../../../source/ics/FunctionIC.html"}>>>
function<<<{"description": "The initial condition function."}>>> = 1-x
variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = sat
[]
[]
[Kernels<<<{"href": "../../../../syntax/Kernels/index.html"}>>>]
[mass0]
type = PorousFlowMassTimeDerivative<<<{"description": "Derivative of fluid-component mass with respect to time. Mass lumping to the nodes is used.", "href": "../../../../source/kernels/PorousFlowMassTimeDerivative.html"}>>>
fluid_component<<<{"description": "The index corresponding to the component for this kernel"}>>> = 0
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pp
[]
[mass1]
type = PorousFlowMassTimeDerivative<<<{"description": "Derivative of fluid-component mass with respect to time. Mass lumping to the nodes is used.", "href": "../../../../source/kernels/PorousFlowMassTimeDerivative.html"}>>>
fluid_component<<<{"description": "The index corresponding to the component for this kernel"}>>> = 1
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = sat
[]
[]
[UserObjects<<<{"href": "../../../../syntax/UserObjects/index.html"}>>>]
[dictator]
type = PorousFlowDictator<<<{"description": "Holds information on the PorousFlow variable names", "href": "../../../../source/userobjects/PorousFlowDictator.html"}>>>
porous_flow_vars<<<{"description": "List of primary variables that are used in the PorousFlow simulation. Jacobian entries involving derivatives wrt these variables will be computed. In single-phase models you will just have one (eg 'pressure'), in two-phase models you will have two (eg 'p_water p_gas', or 'p_water s_water'), etc."}>>> = 'pp sat'
number_fluid_phases<<<{"description": "The number of fluid phases in the simulation"}>>> = 2
number_fluid_components<<<{"description": "The number of fluid components in the simulation"}>>> = 2
[]
[pc]
type = PorousFlowCapillaryPressureConst<<<{"description": "Constant capillary pressure", "href": "../../../../source/userobjects/PorousFlowCapillaryPressureConst.html"}>>>
pc<<<{"description": "Constant capillary pressure (Pa). Default is 0"}>>> = 0
[]
[]
[FluidProperties<<<{"href": "../../../../syntax/FluidProperties/index.html"}>>>]
[simple_fluid0]
type = SimpleFluidProperties<<<{"description": "Fluid properties for a simple fluid with a constant bulk density", "href": "../../../../source/fluidproperties/SimpleFluidProperties.html"}>>>
bulk_modulus<<<{"description": "Constant bulk modulus (Pa)"}>>> = 1
density0<<<{"description": "Density at zero pressure and zero temperature"}>>> = 1
thermal_expansion<<<{"description": "Constant coefficient of thermal expansion (1/K)"}>>> = 0
[]
[simple_fluid1]
type = SimpleFluidProperties<<<{"description": "Fluid properties for a simple fluid with a constant bulk density", "href": "../../../../source/fluidproperties/SimpleFluidProperties.html"}>>>
bulk_modulus<<<{"description": "Constant bulk modulus (Pa)"}>>> = 1
density0<<<{"description": "Density at zero pressure and zero temperature"}>>> = 0.1
thermal_expansion<<<{"description": "Constant coefficient of thermal expansion (1/K)"}>>> = 0
[]
[]
[Materials<<<{"href": "../../../../syntax/Materials/index.html"}>>>]
[temperature]
type = PorousFlowTemperature<<<{"description": "Material to provide temperature at the quadpoints or nodes and derivatives of it with respect to the PorousFlow variables", "href": "../../../../source/materials/PorousFlowTemperature.html"}>>>
[]
[ppss]
type = PorousFlow2PhasePS<<<{"description": "This Material calculates the 2 porepressures and the 2 saturations in a 2-phase situation, and derivatives of these with respect to the PorousFlowVariables.", "href": "../../../../source/materials/PorousFlow2PhasePS.html"}>>>
phase0_porepressure<<<{"description": "Variable that is the porepressure of phase 0 (the liquid phase)"}>>> = pp
phase1_saturation<<<{"description": "Variable that is the saturation of phase 1 (the gas phase)"}>>> = sat
capillary_pressure<<<{"description": "Name of the UserObject defining the capillary pressure"}>>> = pc
[]
[massfrac]
type = PorousFlowMassFraction<<<{"description": "This Material forms a std::vector<std::vector ...> of mass-fractions out of the individual mass fractions", "href": "../../../../source/materials/PorousFlowMassFraction.html"}>>>
mass_fraction_vars<<<{"description": "List of variables that represent the mass fractions. Format is 'f_ph0^c0 f_ph0^c1 f_ph0^c2 ... f_ph0^c(N-2) f_ph1^c0 f_ph1^c1 fph1^c2 ... fph1^c(N-2) ... fphP^c0 f_phP^c1 fphP^c2 ... fphP^c(N-2)' where N=num_components and P=num_phases, and it is assumed that f_ph^c(N-1)=1-sum(f_ph^c,{c,0,N-2}) so that f_ph^c(N-1) need not be given. If no variables are provided then num_phases=1=num_components."}>>> = 'massfrac_ph0_sp0 massfrac_ph1_sp0'
[]
[simple_fluid0]
type = PorousFlowSingleComponentFluid<<<{"description": "This Material calculates fluid properties at the quadpoints or nodes for a single component fluid", "href": "../../../../source/materials/PorousFlowSingleComponentFluid.html"}>>>
fp<<<{"description": "The name of the user object for fluid properties"}>>> = simple_fluid0
phase<<<{"description": "The phase number"}>>> = 0
[]
[simple_fluid1]
type = PorousFlowSingleComponentFluid<<<{"description": "This Material calculates fluid properties at the quadpoints or nodes for a single component fluid", "href": "../../../../source/materials/PorousFlowSingleComponentFluid.html"}>>>
fp<<<{"description": "The name of the user object for fluid properties"}>>> = simple_fluid1
phase<<<{"description": "The phase number"}>>> = 1
[]
[porosity]
type = PorousFlowPorosityConst<<<{"description": "This Material calculates the porosity assuming it is constant", "href": "../../../../source/materials/PorousFlowPorosityConst.html"}>>>
porosity<<<{"description": "The porosity (assumed indepenent of porepressure, temperature, strain, etc, for this material). This should be a real number, or a constant monomial variable (not a linear lagrange or other kind of variable)."}>>> = 0.1
[]
[]
[Postprocessors<<<{"href": "../../../../syntax/Postprocessors/index.html"}>>>]
[comp0_phase0_mass]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 0
phase<<<{"description": "The index of the fluid phase that this Postprocessor is restricted to. Multiple indices can be entered"}>>> = 0
[]
[comp0_phase1_mass]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 0
phase<<<{"description": "The index of the fluid phase that this Postprocessor is restricted to. Multiple indices can be entered"}>>> = 1
[]
[comp0_total_mass]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 0
[]
[comp1_phase0_mass]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 1
phase<<<{"description": "The index of the fluid phase that this Postprocessor is restricted to. Multiple indices can be entered"}>>> = 0
[]
[comp1_phase1_mass]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 1
phase<<<{"description": "The index of the fluid phase that this Postprocessor is restricted to. Multiple indices can be entered"}>>> = 1
[]
[comp1_total_mass]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 1
[]
[comp1_phase1_threshold_mass]
type = PorousFlowFluidMass<<<{"description": "Calculates the mass of a fluid component in a region", "href": "../../../../source/postprocessors/PorousFlowFluidMass.html"}>>>
fluid_component<<<{"description": "The index corresponding to the fluid component that this Postprocessor acts on"}>>> = 1
phase<<<{"description": "The index of the fluid phase that this Postprocessor is restricted to. Multiple indices can be entered"}>>> = 1
saturation_threshold<<<{"description": "The saturation threshold below which the mass is calculated for a specific phase. Default is 1.0. Note: only one phase_index can be entered"}>>> = 0.6
[]
[]
[Executioner<<<{"href": "../../../../syntax/Executioner/index.html"}>>>]
type = Transient
solve_type = Newton
nl_abs_tol = 1e-16
dt = 1
end_time = 1
[]
[Outputs<<<{"href": "../../../../syntax/Outputs/index.html"}>>>]
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'timestep_end'
file_base<<<{"description": "Common file base name to be utilized with all output objects"}>>> = mass06
csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
[]
(moose/modules/porous_flow/test/tests/mass_conservation/mass06.i)