- boundaryThe list of boundary IDs from the mesh where this boundary condition applies
C++ Type:std::vector<BoundaryName>
Description:The list of boundary IDs from the mesh where this boundary condition applies
AreaPostprocessor
The AreaPostprocesor is a SidePostprocessor that computes the area or dimension - 1 "volume" of a given side of the mesh. This postprocessor may be applied to one or more boundaries simulatenously, the latter case produces a total area is output.
Description and Syntax
Computes the "area" or dimension - 1 "volume" of a given boundary or boundaries in your mesh.
Input Parameters
- execute_onTIMESTEP_ENDThe list of flag(s) indicating when this object should be executed, the available options include NONE, INITIAL, LINEAR, NONLINEAR, TIMESTEP_END, TIMESTEP_BEGIN, FINAL, CUSTOM.
Default:TIMESTEP_END
C++ Type:ExecFlagEnum
Options:NONE, INITIAL, LINEAR, NONLINEAR, TIMESTEP_END, TIMESTEP_BEGIN, FINAL, CUSTOM, TRANSFER
Description:The list of flag(s) indicating when this object should be executed, the available options include NONE, INITIAL, LINEAR, NONLINEAR, TIMESTEP_END, TIMESTEP_BEGIN, FINAL, CUSTOM.
Optional Parameters
- allow_duplicate_execution_on_initialFalseIn the case where this UserObject is depended upon by an initial condition, allow it to be executed twice during the initial setup (once before the IC and again after mesh adaptivity (if applicable).
Default:False
C++ Type:bool
Options:
Description:In the case where this UserObject is depended upon by an initial condition, allow it to be executed twice during the initial setup (once before the IC and again after mesh adaptivity (if applicable).
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Options:
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Options:
Description:Set the enabled status of the MooseObject.
- outputsVector of output names were you would like to restrict the output of variables(s) associated with this object
C++ Type:std::vector<OutputName>
Options:
Description:Vector of output names were you would like to restrict the output of variables(s) associated with this object
- use_displaced_meshFalseWhether 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.
Default:False
C++ Type:bool
Options:
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.
Advanced Parameters
Input Files
- (modules/porous_flow/examples/restart/gas_injection_new_mesh.i)
- (test/tests/postprocessors/area_pp/area_pp.i)
- (test/tests/postprocessors/geometry/3d_geometry.i)
- (modules/porous_flow/test/tests/fluidstate/coldwater_injection_radial.i)
- (modules/richards/test/tests/sinks/s04.i)
- (test/tests/postprocessors/geometry/2d_geometry.i)
- (modules/richards/test/tests/sinks/s_fu_04.i)
- (modules/porous_flow/examples/restart/gas_injection.i)
References
(modules/porous_flow/examples/restart/gas_injection_new_mesh.i)
# Using the results from the equilibrium run to provide the initial condition for
# porepressure, we now inject a gas phase into the brine-saturated reservoir. In this
# example, the mesh is not identical to the mesh used in gravityeq.i. Rather, it is
# generated so that it is more refined near the injection boundary and at the top of
# the model, as that is where the gas plume will be present.
#
# To use the hydrostatic pressure calculated using the gravity equilibrium run as the initial
# condition for the pressure, a SolutionUserObject is used, along with a SolutionFunction to
# interpolate the pressure from the gravity equilibrium run to the initial condition for liqiud
# porepressure in this example.
#
# Even though the gravity equilibrium is established using a 2D mesh, in this example,
# we use a mesh shifted 0.1 m to the right and rotate it about the Y axis to make a 2D radial
# model.
#
# Methane injection takes place over the surface of the hole created by rotating the mesh,
# and hence the injection area is 2 pi r h. We can calculate this using an AreaPostprocessor,
# and then use this in a ParsedFunction to calculate the injection rate so that 10 kg/s of
# methane is injected.
#
# Note: as this example uses the results from a previous simulation, gravityeq.i MUST be
# run before running this input file.
[Mesh]
type = GeneratedMesh
dim = 2
ny = 25
nx = 50
ymax = 100
xmin = 0.1
xmax = 5000
bias_x = 1.05
bias_y = 0.95
[]
[Problem]
coord_type = RZ
rz_coord_axis = Y
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 -9.81 0'
temperature_unit = Celsius
[]
[Variables]
[pp_liq]
[]
[sat_gas]
initial_condition = 0
[]
[]
[ICs]
[ppliq_ic]
type = FunctionIC
variable = pp_liq
function = ppliq_ic
[]
[]
[AuxVariables]
[temperature]
initial_condition = 50
[]
[xnacl]
initial_condition = 0.1
[]
[brine_density]
family = MONOMIAL
order = CONSTANT
[]
[methane_density]
family = MONOMIAL
order = CONSTANT
[]
[massfrac_ph0_sp0]
initial_condition = 1
[]
[massfrac_ph1_sp0]
initial_condition = 0
[]
[pp_gas]
family = MONOMIAL
order = CONSTANT
[]
[sat_liq]
family = MONOMIAL
order = CONSTANT
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
variable = pp_liq
[]
[flux0]
type = PorousFlowAdvectiveFlux
variable = pp_liq
[]
[mass1]
type = PorousFlowMassTimeDerivative
variable = sat_gas
fluid_component = 1
[]
[flux1]
type = PorousFlowAdvectiveFlux
variable = sat_gas
fluid_component = 1
[]
[]
[AuxKernels]
[brine_density]
type = PorousFlowPropertyAux
property = density
variable = brine_density
execute_on = 'initial timestep_end'
[]
[methane_density]
type = PorousFlowPropertyAux
property = density
variable = methane_density
phase = 1
execute_on = 'initial timestep_end'
[]
[pp_gas]
type = PorousFlowPropertyAux
property = pressure
phase = 1
variable = pp_gas
execute_on = 'initial timestep_end'
[]
[sat_liq]
type = PorousFlowPropertyAux
property = saturation
variable = sat_liq
execute_on = 'initial timestep_end'
[]
[]
[BCs]
[gas_injection]
type = PorousFlowSink
boundary = left
variable = sat_gas
flux_function = injection_rate
fluid_phase = 1
[]
[brine_out]
type = PorousFlowPiecewiseLinearSink
boundary = right
variable = pp_liq
multipliers = '0 1e9'
pt_vals = '0 1e9'
fluid_phase = 0
flux_function = 1e-6
use_mobility = true
use_relperm = true
mass_fraction_component = 0
[]
[]
[Functions]
[injection_rate]
type = ParsedFunction
vals = injection_area
vars = area
value = '-1/area'
[]
[ppliq_ic]
type = SolutionFunction
solution = soln
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp_liq sat_gas'
number_fluid_phases = 2
number_fluid_components = 2
[]
[pc]
type = PorousFlowCapillaryPressureVG
alpha = 1e-5
m = 0.5
sat_lr = 0.2
pc_max = 1e7
[]
[soln]
type = SolutionUserObject
mesh = gravityeq_out.e
system_variables = porepressure
[]
[]
[Modules]
[FluidProperties]
[brine]
type = BrineFluidProperties
[]
[methane]
type = MethaneFluidProperties
[]
[methane_tab]
type = TabulatedFluidProperties
fp = methane
save_file = false
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
temperature = temperature
[]
[ps]
type = PorousFlow2PhasePS
phase0_porepressure = pp_liq
phase1_saturation = sat_gas
capillary_pressure = pc
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = 'massfrac_ph0_sp0 massfrac_ph1_sp0'
[]
[brine]
type = PorousFlowBrine
compute_enthalpy = false
compute_internal_energy = false
xnacl = xnacl
phase = 0
[]
[methane]
type = PorousFlowSingleComponentFluid
compute_enthalpy = false
compute_internal_energy = false
fp = methane_tab
phase = 1
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-13 0 0 0 5e-14 0 0 0 1e-13'
[]
[relperm_liq]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
s_res = 0.2
sum_s_res = 0.3
[]
[relperm_gas]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 1
s_res = 0.1
sum_s_res = 0.3
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type'
petsc_options_value = ' asm lu NONZERO'
[]
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 1e8
nl_abs_tol = 1e-12
nl_rel_tol = 1e-06
nl_max_its = 20
dtmax = 1e6
[TimeStepper]
type = IterationAdaptiveDT
dt = 1e1
growth_factor = 1.5
[]
[]
[Postprocessors]
[mass_ph0]
type = PorousFlowFluidMass
fluid_component = 0
execute_on = 'initial timestep_end'
[]
[mass_ph1]
type = PorousFlowFluidMass
fluid_component = 1
execute_on = 'initial timestep_end'
[]
[injection_area]
type = AreaPostprocessor
boundary = left
execute_on = initial
[]
[]
[Outputs]
execute_on = 'initial timestep_end'
exodus = true
perf_graph = true
[]
(test/tests/postprocessors/area_pp/area_pp.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
ny = 4
xmax = 1.2
ymax = 2.3
[]
[Variables]
[./u]
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = left
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = right
value = 1
[../]
[]
[Postprocessors]
[./right]
type = AreaPostprocessor
boundary = 'right'
execute_on = 'initial timestep_end'
[../]
[./bottom]
type = AreaPostprocessor
boundary = 'bottom'
execute_on = 'initial timestep_end'
[../]
[./all]
type = AreaPostprocessor
boundary = 'left right bottom top'
execute_on = 'initial timestep_end'
[../]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
csv = true
[]
(test/tests/postprocessors/geometry/3d_geometry.i)
radius = 0.5
inner_box_length = 2.2
outer_box_length = 3
depth = 0.4
sides = 28
alpha = ${fparse 2 * pi / ${sides}}
perimeter_correction = ${fparse ${alpha} / 2 / sin(alpha / 2)}
area_correction = ${fparse alpha / sin(alpha)}
[Mesh]
file = 3d.e
construct_side_list_from_node_list = true
[]
[Problem]
solve = false
[]
[Executioner]
type = Steady
[]
[Postprocessors]
[./circle_side_area]
type = AreaPostprocessor
boundary = circle_side
[../]
[./inside_side_area]
type = AreaPostprocessor
boundary = inside_side
[../]
[./outside_side_area]
type = AreaPostprocessor
boundary = outside_side
[../]
[./circle_volume]
type = VolumePostprocessor
block = circle
[../]
[./inside_volume]
type = VolumePostprocessor
block = inside
[../]
[./outside_volume]
type = VolumePostprocessor
block = outside
[../]
[./total_volume]
type = VolumePostprocessor
block = 'circle inside outside'
[../]
[./circle_side_area_exact]
type = FunctionValuePostprocessor
function = 'circle_side_area_exact'
[../]
[./inside_side_area_exact]
type = FunctionValuePostprocessor
function = 'inside_side_area_exact'
[../]
[./outside_side_area_exact]
type = FunctionValuePostprocessor
function = 'outside_side_area_exact'
[../]
[./circle_volume_exact]
type = FunctionValuePostprocessor
function = 'circle_volume_exact'
[../]
[./inside_volume_exact]
type = FunctionValuePostprocessor
function = 'inside_volume_exact'
[../]
[./outside_volume_exact]
type = FunctionValuePostprocessor
function = 'outside_volume_exact'
[../]
[./total_volume_exact]
type = FunctionValuePostprocessor
function = 'total_volume_exact'
[../]
[]
[Functions]
[./circle_side_area_exact]
type = ParsedFunction
value = '2 * pi * ${radius} / ${perimeter_correction} * ${depth}'
[../]
[./inside_side_area_exact]
type = ParsedFunction
value = '${inner_box_length} * ${depth} * 4'
[../]
[./outside_side_area_exact]
type = ParsedFunction
value = '${outer_box_length} * ${depth} * 4'
[../]
[./circle_volume_exact]
type = ParsedFunction
value = 'pi * ${radius}^2 * ${depth} / ${area_correction}'
[../]
[./inside_volume_exact]
type = ParsedFunction
value = '${inner_box_length}^2 * ${depth} - pi * ${radius}^2 * ${depth} / ${area_correction}'
[../]
[./outside_volume_exact]
type = ParsedFunction
value = '${outer_box_length}^2 * ${depth} - ${inner_box_length}^2 * ${depth}'
[../]
[./total_volume_exact]
type = ParsedFunction
value = '${outer_box_length}^2 * ${depth}'
[../]
[]
[Outputs]
csv = true
[]
(modules/porous_flow/test/tests/fluidstate/coldwater_injection_radial.i)
# Cold water injection into 1D radial hot reservoir (Avdonin, 1964)
#
# To generate results presented in documentation for this problem,
# set xmax = 1000 and nx = 200 in the Mesh block, and dtmax = 1e4
# and end_time = 1e6 in the Executioner block.
[Mesh]
type = GeneratedMesh
dim = 1
nx = 50
xmin = 0.1
xmax = 5
bias_x = 1.05
[]
[Problem]
rz_coord_axis = Y
coord_type = RZ
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[AuxVariables]
[temperature]
order = CONSTANT
family = MONOMIAL
[]
[]
[AuxKernels]
[temperature]
type = PorousFlowPropertyAux
variable = temperature
property = temperature
execute_on = 'initial timestep_end'
[]
[]
[Variables]
[pliquid]
initial_condition = 5e6
[]
[h]
scaling = 1e-6
[]
[]
[ICs]
[hic]
type = PorousFlowFluidPropertyIC
variable = h
porepressure = pliquid
property = enthalpy
temperature = 170
temperature_unit = Celsius
fp = water
[]
[]
[Functions]
[injection_rate]
type = ParsedFunction
vals = injection_area
vars = area
value = '-0.1/area'
[]
[]
[BCs]
[source]
type = PorousFlowSink
variable = pliquid
flux_function = injection_rate
boundary = left
[]
[pright]
type = DirichletBC
variable = pliquid
value = 5e6
boundary = right
[]
[hleft]
type = DirichletBC
variable = h
value = 678.52e3
boundary = left
[]
[hright]
type = DirichletBC
variable = h
value = 721.4e3
boundary = right
[]
[]
[Kernels]
[mass]
type = PorousFlowMassTimeDerivative
variable = pliquid
[]
[massflux]
type = PorousFlowAdvectiveFlux
variable = pliquid
[]
[heat]
type = PorousFlowEnergyTimeDerivative
variable = h
[]
[heatflux]
type = PorousFlowHeatAdvection
variable = h
[]
[heatcond]
type = PorousFlowHeatConduction
variable = h
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pliquid h'
number_fluid_phases = 2
number_fluid_components = 1
[]
[pc]
type = PorousFlowCapillaryPressureVG
pc_max = 1e6
sat_lr = 0.1
m = 0.5
alpha = 1e-5
[]
[fs]
type = PorousFlowWaterVapor
water_fp = water
capillary_pressure = pc
[]
[]
[Modules]
[FluidProperties]
[water]
type = Water97FluidProperties
[]
[]
[]
[Materials]
[watervapor]
type = PorousFlowFluidStateSingleComponent
porepressure = pliquid
enthalpy = h
temperature_unit = Celsius
capillary_pressure = pc
fluid_state = fs
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.2
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1.8e-11 0 0 0 1.8e-11 0 0 0 1.8e-11'
[]
[relperm_water]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
s_res = 0.1
sum_s_res = 0.1
[]
[relperm_gas]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 1
sum_s_res = 0.1
[]
[internal_energy]
type = PorousFlowMatrixInternalEnergy
density = 2900
specific_heat_capacity = 740
[]
[rock_thermal_conductivity]
type = PorousFlowThermalConductivityIdeal
dry_thermal_conductivity = '20 0 0 0 20 0 0 0 20'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
end_time = 1e3
nl_abs_tol = 1e-8
[TimeStepper]
type = IterationAdaptiveDT
dt = 100
[]
[]
[Postprocessors]
[injection_area]
type = AreaPostprocessor
boundary = left
execute_on = initial
[]
[]
[VectorPostprocessors]
[line]
type = ElementValueSampler
sort_by = x
variable = temperature
execute_on = 'initial timestep_end'
[]
[]
[Outputs]
perf_graph = true
[csv]
type = CSV
execute_on = final
[]
[]
(modules/richards/test/tests/sinks/s04.i)
# apply a total flux (in kg/s) to two boundaries
# and check that it removes the correct amount of fluid
[Mesh]
type = GeneratedMesh
dim = 2
nx = 1
ny = 3
xmin = 0
xmax = 1
ymin = 0
ymax = 4
[]
[GlobalParams]
richardsVarNames_UO = PPNames
density_UO = DensityConstBulk
relperm_UO = RelPermPower
SUPG_UO = SUPGstandard
sat_UO = Saturation
seff_UO = SeffVG
viscosity = 1E-3
gravity = '-1 0 0'
[]
[UserObjects]
[./PPNames]
type = RichardsVarNames
richards_vars = pressure
[../]
[./DensityConstBulk]
type = RichardsDensityConstBulk
dens0 = 1
bulk_mod = 1
[../]
[./SeffVG]
type = RichardsSeff1VG
m = 0.5
al = 1
[../]
[./RelPermPower]
type = RichardsRelPermPower
simm = 0.0
n = 2
[../]
[./Saturation]
type = RichardsSat
s_res = 0.1
sum_s_res = 0.2
[../]
[./SUPGstandard]
type = RichardsSUPGstandard
p_SUPG = 0.1
[../]
[]
[Variables]
[./pressure]
[../]
[]
[ICs]
[./pressure]
type = ConstantIC
variable = pressure
value = 2
[../]
[]
[Postprocessors]
[./area_left]
type = AreaPostprocessor
boundary = left
execute_on = initial
[../]
[./area_right]
type = AreaPostprocessor
boundary = right
execute_on = initial
[../]
[./mass_fin]
type = RichardsMass
variable = pressure
execute_on = 'initial timestep_end'
[../]
[./p0]
type = PointValue
point = '0 0 0'
variable = pressure
execute_on = 'initial timestep_end'
[../]
[]
[BCs]
[./left_flux]
type = RichardsPiecewiseLinearSink
boundary = left
pressures = '0'
bare_fluxes = '0.1'
variable = pressure
use_mobility = false
use_relperm = false
area_pp = area_left
[../]
[./right_flux]
type = RichardsPiecewiseLinearSink
boundary = right
pressures = '0'
bare_fluxes = '0.1'
variable = pressure
use_mobility = false
use_relperm = false
area_pp = area_right
[../]
[]
[Kernels]
active = 'richardst'
[./richardst]
type = RichardsMassChange
variable = pressure
[../]
[]
[Materials]
[./rock]
type = RichardsMaterial
block = 0
mat_porosity = 0.1
mat_permeability = '1E-5 0 0 0 1E-5 0 0 0 1E-5'
linear_shape_fcns = true
[../]
[]
[Preconditioning]
[./andy]
type = SMP
full = true
petsc_options_iname = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it'
petsc_options_value = 'bcgs bjacobi 1E-12 1E-10 10000'
[../]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 1
end_time = 13
[]
[Outputs]
file_base = s04
csv = true
[]
(test/tests/postprocessors/geometry/2d_geometry.i)
radius = 0.5
inner_box_length = 2.2
outer_box_length = 3
sides = 16
alpha = ${fparse 2 * pi / ${sides}}
perimeter_correction = ${fparse alpha / 2 / sin(alpha / 2)}
area_correction = ${fparse alpha / sin(alpha)}
[Mesh]
file = 2d.e
construct_side_list_from_node_list = true
[]
[Variables]
[./u]
initial_condition = 1
block = circle
[../]
[./v]
initial_condition = 2
block = 'inside outside'
[../]
[]
[Kernels]
[./diff_u]
type = Diffusion
variable = u
[../]
[./diff_v]
type = Diffusion
variable = v
[../]
[]
[BCs]
[./circle]
type = DirichletBC
variable = u
boundary = circle_side_wrt_inside
value = 2
[../]
[./inner]
type = DirichletBC
variable = v
boundary = circle_side_wrt_circle
value = 4
[../]
[./outer]
type = DirichletBC
variable = v
boundary = outside_side
value = 6
[../]
[]
[Executioner]
type = Steady
[]
[Postprocessors]
[./u_avg]
type = ElementAverageValue
variable = u
block = circle
[../]
[./v_avg]
type = ElementAverageValue
variable = v
block = 'inside outside'
[../]
[./circle_perimeter_wrt_circle]
type = AreaPostprocessor
boundary = circle_side_wrt_circle
[../]
[./circle_perimeter_wrt_inside]
type = AreaPostprocessor
boundary = circle_side_wrt_inside
[../]
[./inside_perimeter_wrt_inside]
type = AreaPostprocessor
boundary = inside_side_wrt_inside
[../]
[./inside_perimeter_wrt_outside]
type = AreaPostprocessor
boundary = inside_side_wrt_outside
[../]
[./outside_perimeter]
type = AreaPostprocessor
boundary = outside_side
[../]
[./circle_area]
type = VolumePostprocessor
block = circle
[../]
[./inside_area]
type = VolumePostprocessor
block = inside
[../]
[./outside_area]
type = VolumePostprocessor
block = outside
[../]
[./total_area]
type = VolumePostprocessor
block = 'circle inside outside'
[../]
[./circle_perimeter_exact]
type = FunctionValuePostprocessor
function = 'circle_perimeter_exact'
[../]
[./inside_perimeter_exact]
type = FunctionValuePostprocessor
function = 'inside_perimeter_exact'
[../]
[./outside_perimeter_exact]
type = FunctionValuePostprocessor
function = 'outside_perimeter_exact'
[../]
[./circle_area_exact]
type = FunctionValuePostprocessor
function = 'circle_area_exact'
[../]
[./inside_area_exact]
type = FunctionValuePostprocessor
function = 'inside_area_exact'
[../]
[./outside_area_exact]
type = FunctionValuePostprocessor
function = 'outside_area_exact'
[../]
[./total_area_exact]
type = FunctionValuePostprocessor
function = 'total_area_exact'
[../]
[]
[Functions]
[./circle_perimeter_exact]
type = ParsedFunction
value = '2 * pi * ${radius} / ${perimeter_correction}'
[../]
[./inside_perimeter_exact]
type = ParsedFunction
value = '${inner_box_length} * 4'
[../]
[./outside_perimeter_exact]
type = ParsedFunction
value = '${outer_box_length} * 4'
[../]
[./circle_area_exact]
type = ParsedFunction
value = 'pi * ${radius}^2 / ${area_correction}'
[../]
[./inside_area_exact]
type = ParsedFunction
value = '${inner_box_length}^2 - pi * ${radius}^2 / ${area_correction}'
[../]
[./outside_area_exact]
type = ParsedFunction
value = '${outer_box_length}^2 - ${inner_box_length}^2'
[../]
[./total_area_exact]
type = ParsedFunction
value = '${outer_box_length}^2'
[../]
[]
[Outputs]
csv = true
[]
(modules/richards/test/tests/sinks/s_fu_04.i)
# apply a total flux (in kg/s) to two boundaries
# and check that it removes the correct amount of fluid
# fully-upwind sink
[Mesh]
type = GeneratedMesh
dim = 2
nx = 1
ny = 3
xmin = 0
xmax = 1
ymin = 0
ymax = 4
[]
[GlobalParams]
richardsVarNames_UO = PPNames
density_UO = DensityConstBulk
relperm_UO = RelPermPower
SUPG_UO = SUPGstandard
sat_UO = Saturation
seff_UO = SeffVG
viscosity = 1E-3
gravity = '-1 0 0'
[]
[UserObjects]
[./PPNames]
type = RichardsVarNames
richards_vars = pressure
[../]
[./DensityConstBulk]
type = RichardsDensityConstBulk
dens0 = 1
bulk_mod = 1
[../]
[./SeffVG]
type = RichardsSeff1VG
m = 0.5
al = 1 # same deal with PETSc constant state
[../]
[./RelPermPower]
type = RichardsRelPermPower
simm = 0.0
n = 2
[../]
[./Saturation]
type = RichardsSat
s_res = 0.1
sum_s_res = 0.2
[../]
[./SUPGstandard]
type = RichardsSUPGstandard
p_SUPG = 0.1
[../]
[]
[Variables]
[./pressure]
[../]
[]
[ICs]
[./pressure]
type = ConstantIC
variable = pressure
value = 2
[../]
[]
[Postprocessors]
[./area_left]
type = AreaPostprocessor
boundary = left
execute_on = initial
[../]
[./area_right]
type = AreaPostprocessor
boundary = right
execute_on = initial
[../]
[./mass_fin]
type = RichardsMass
variable = pressure
execute_on = 'initial timestep_end'
[../]
[./p0]
type = PointValue
point = '0 0 0'
variable = pressure
execute_on = 'initial timestep_end'
[../]
[]
[BCs]
[./left_flux]
type = RichardsPiecewiseLinearSink
boundary = left
pressures = '0'
bare_fluxes = '0.1'
variable = pressure
use_mobility = false
use_relperm = false
area_pp = area_left
fully_upwind = true
[../]
[./right_flux]
type = RichardsPiecewiseLinearSink
boundary = right
pressures = '0'
bare_fluxes = '0.1'
variable = pressure
use_mobility = false
use_relperm = false
area_pp = area_right
fully_upwind = true
[../]
[]
[Kernels]
active = 'richardst'
[./richardst]
type = RichardsMassChange
variable = pressure
[../]
[]
[Materials]
[./rock]
type = RichardsMaterial
block = 0
mat_porosity = 0.1
mat_permeability = '1E-5 0 0 0 1E-5 0 0 0 1E-5'
linear_shape_fcns = true
[../]
[]
[Preconditioning]
[./andy]
type = SMP
full = true
petsc_options_iname = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it'
petsc_options_value = 'bcgs bjacobi 1E-12 1E-10 10000'
[../]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 1
end_time = 13
[]
[Outputs]
file_base = s_fu_04
csv = true
[]
(modules/porous_flow/examples/restart/gas_injection.i)
# Using the results from the equilibrium run to provide the initial condition for
# porepressure, we now inject a gas phase into the brine-saturated reservoir. In this
# example, where the mesh used is identical to the mesh used in gravityeq.i, we can use
# the basic restart capability by simply setting the initial condition for porepressure
# using the results from gravityeq.i.
#
# Even though the gravity equilibrium is established using a 2D mesh, in this example,
# we shift the mesh 0.1 m to the right and rotate it about the Y axis to make a 2D radial
# model.
#
# Methane injection takes place over the surface of the hole created by rotating the mesh,
# and hence the injection area is 2 pi r h. We can calculate this using an AreaPostprocessor,
# and then use this in a ParsedFunction to calculate the injection rate so that 10 kg/s of
# methane is injected.
#
# Results can be improved by uniformly refining the initial mesh.
#
# Note: as this example uses the results from a previous simulation, gravityeq.i MUST be
# run before running this input file.
[Mesh]
uniform_refine = 1
[file]
type = FileMeshGenerator
file = gravityeq_out.e
[]
[translate]
type = TransformGenerator
transform = TRANSLATE
vector_value = '0.1 0 0'
input = file
[]
[]
[Problem]
coord_type = RZ
rz_coord_axis = Y
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 -9.81 0'
temperature_unit = Celsius
[]
[Variables]
[pp_liq]
initial_from_file_var = porepressure
[]
[sat_gas]
initial_condition = 0
[]
[]
[AuxVariables]
[temperature]
initial_condition = 50
[]
[xnacl]
initial_condition = 0.1
[]
[brine_density]
family = MONOMIAL
order = CONSTANT
[]
[methane_density]
family = MONOMIAL
order = CONSTANT
[]
[massfrac_ph0_sp0]
initial_condition = 1
[]
[massfrac_ph1_sp0]
initial_condition = 0
[]
[pp_gas]
family = MONOMIAL
order = CONSTANT
[]
[sat_liq]
family = MONOMIAL
order = CONSTANT
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
variable = pp_liq
[]
[flux0]
type = PorousFlowAdvectiveFlux
variable = pp_liq
[]
[mass1]
type = PorousFlowMassTimeDerivative
variable = sat_gas
fluid_component = 1
[]
[flux1]
type = PorousFlowAdvectiveFlux
variable = sat_gas
fluid_component = 1
[]
[]
[AuxKernels]
[brine_density]
type = PorousFlowPropertyAux
property = density
variable = brine_density
execute_on = 'initial timestep_end'
[]
[methane_density]
type = PorousFlowPropertyAux
property = density
variable = methane_density
phase = 1
execute_on = 'initial timestep_end'
[]
[pp_gas]
type = PorousFlowPropertyAux
property = pressure
phase = 1
variable = pp_gas
execute_on = 'initial timestep_end'
[]
[sat_liq]
type = PorousFlowPropertyAux
property = saturation
variable = sat_liq
execute_on = 'initial timestep_end'
[]
[]
[BCs]
[gas_injection]
type = PorousFlowSink
boundary = left
variable = sat_gas
flux_function = injection_rate
fluid_phase = 1
[]
[brine_out]
type = PorousFlowPiecewiseLinearSink
boundary = right
variable = pp_liq
multipliers = '0 1e9'
pt_vals = '0 1e9'
fluid_phase = 0
flux_function = 1e-6
use_mobility = true
[]
[]
[Functions]
[injection_rate]
type = ParsedFunction
vals = injection_area
vars = area
value = '-10/area'
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp_liq sat_gas'
number_fluid_phases = 2
number_fluid_components = 2
[]
[pc]
type = PorousFlowCapillaryPressureVG
alpha = 1e-5
m = 0.5
sat_lr = 0.2
[]
[]
[Modules]
[FluidProperties]
[brine]
type = BrineFluidProperties
[]
[methane]
type = MethaneFluidProperties
[]
[methane_tab]
type = TabulatedFluidProperties
fp = methane
save_file = false
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
temperature = temperature
[]
[ps]
type = PorousFlow2PhasePS
phase0_porepressure = pp_liq
phase1_saturation = sat_gas
capillary_pressure = pc
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = 'massfrac_ph0_sp0 massfrac_ph1_sp0'
[]
[brine]
type = PorousFlowBrine
compute_enthalpy = false
compute_internal_energy = false
xnacl = xnacl
phase = 0
[]
[methane]
type = PorousFlowSingleComponentFluid
compute_enthalpy = false
compute_internal_energy = false
fp = methane_tab
phase = 1
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-13 0 0 0 1e-13 0 0 0 1e-13'
[]
[relperm_liq]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
s_res = 0.2
sum_s_res = 0.3
[]
[relperm_gas]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 1
s_res = 0.1
sum_s_res = 0.3
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type'
petsc_options_value = ' asm lu NONZERO'
[]
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 1e8
nl_abs_tol = 1e-12
nl_rel_tol = 1e-06
nl_max_its = 20
dtmax = 1e6
[TimeStepper]
type = IterationAdaptiveDT
dt = 1e1
[]
[]
[Postprocessors]
[mass_ph0]
type = PorousFlowFluidMass
fluid_component = 0
execute_on = 'initial timestep_end'
[]
[mass_ph1]
type = PorousFlowFluidMass
fluid_component = 1
execute_on = 'initial timestep_end'
[]
[injection_area]
type = AreaPostprocessor
boundary = left
execute_on = initial
[]
[]
[Outputs]
execute_on = 'initial timestep_end'
exodus = true
perf_graph = true
checkpoint = true
[]