LinearFVMomentumBoussinesq
This kernel adds the contributions of the Boussinesq buoyancy treatment for density through a force/source term to the right hand side of the momentum equation system for the finite volume SIMPLE segregated solver SIMPLE. The Boussinesq buoyancy treatment is applicable for low changes in density, and assumes constant density value in all other equation terms.
This term is described by present in the momentum equation conservation when describing an incompressible fluid, where is the reference density, is the thermal expansion coefficient, is the gravity vector, is the temperature, and is a reference temperature. The Boussinesq buoyancy model assumes the changes in density as a function of temperature are linear and relevant only in the buoyant force term of the equation system. The Boussinesq kernel allows for modeling natural convection.
This term deals only with the force due to the variation in density , with the fluid density being . Thus, with no extra added terms to the conventional incompressible Navier Stokes equations, the system will solve for the total pressure minus the hydrostatic pressure. For natural convection simulations, it is advisable to compute relevant dimensionless numbers such as the Rayleigh number or the Richardson number to decide on the need for turbulence models, mesh refinement and stability considerations.
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/2d/2d-boussinesq.i)
mu = 2.6
rho = 1.0
advected_interp_method = 'average'
cp = 300
k = 20
alpha_b = 1e-4
[Mesh]
[mesh]
type = CartesianMeshGenerator
dim = 2
dx = '1.'
dy = '0.2'
ix = '10'
iy = '5'
[]
[]
[Problem]
linear_sys_names = 'u_system v_system pressure_system energy_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = RhieChowMassFlux
u = vel_x
v = vel_y
pressure = pressure
rho = ${rho}
p_diffusion_kernel = p_diffusion
[]
[]
[Variables]
[vel_x]
type = MooseLinearVariableFVReal
initial_condition = 0.5
solver_sys = u_system
[]
[vel_y]
type = MooseLinearVariableFVReal
solver_sys = v_system
initial_condition = 0.0
[]
[pressure]
type = MooseLinearVariableFVReal
solver_sys = pressure_system
initial_condition = 0.2
[]
[T]
type = MooseLinearVariableFVReal
solver_sys = energy_system
initial_condition = 300
[]
[]
[LinearFVKernels]
[u_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_x
advected_interp_method = ${advected_interp_method}
mu = ${mu}
u = vel_x
v = vel_y
momentum_component = 'x'
rhie_chow_user_object = 'rc'
use_nonorthogonal_correction = false
[]
[v_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_y
advected_interp_method = ${advected_interp_method}
mu = ${mu}
u = vel_x
v = vel_y
momentum_component = 'y'
rhie_chow_user_object = 'rc'
use_nonorthogonal_correction = false
[]
[u_pressure]
type = LinearFVMomentumPressure
variable = vel_x
pressure = pressure
momentum_component = 'x'
[]
[v_pressure]
type = LinearFVMomentumPressure
variable = vel_y
pressure = pressure
momentum_component = 'y'
[]
[u_boussinesq]
type = LinearFVMomentumBoussinesq
variable = vel_x
rho = ${rho}
gravity = '0 -9.81 0'
alpha_name = ${alpha_b}
ref_temperature = 300.0
T_fluid = T
momentum_component = 'x'
[]
[v_boussinesq]
type = LinearFVMomentumBoussinesq
variable = vel_y
rho = ${rho}
gravity = '0 -9.81 0'
alpha_name = ${alpha_b}
ref_temperature = 300.0
T_fluid = T
momentum_component = 'y'
[]
[p_diffusion]
type = LinearFVAnisotropicDiffusion
variable = pressure
diffusion_tensor = Ainv
use_nonorthogonal_correction = false
[]
[HbyA_divergence]
type = LinearFVDivergence
variable = pressure
face_flux = HbyA
force_boundary_execution = true
[]
[h_advection]
type = LinearFVEnergyAdvection
variable = T
advected_quantity = temperature
cp = ${cp}
advected_interp_method = ${advected_interp_method}
rhie_chow_user_object = 'rc'
[]
[conduction]
type = LinearFVDiffusion
variable = T
diffusion_coeff = ${k}
use_nonorthogonal_correction = false
[]
[]
[FunctorMaterials]
[constant_functors]
type = GenericFunctorMaterial
prop_names = 'cp alpha_b'
prop_values = '${cp} ${alpha_b}'
[]
[]
[LinearFVBCs]
[inlet-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_x
functor = '1.1'
[]
[inlet-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_y
functor = '0.0'
[]
[walls-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom'
variable = vel_x
functor = 0.0
[]
[walls-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom'
variable = vel_y
functor = 0.0
[]
[outlet_p]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'right'
variable = pressure
functor = 1.4
[]
[outlet_u]
type = LinearFVAdvectionDiffusionOutflowBC
variable = vel_x
use_two_term_expansion = false
boundary = right
[]
[outlet_v]
type = LinearFVAdvectionDiffusionOutflowBC
variable = vel_y
use_two_term_expansion = false
boundary = right
[]
[inlet_top_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T
functor = ${fparse 300.0}
boundary = 'left top'
[]
[bottom_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T
functor = ${fparse 350.0}
boundary = bottom
[]
[outlet_T]
type = LinearFVAdvectionDiffusionOutflowBC
variable = T
use_two_term_expansion = false
boundary = right
[]
[]
[Executioner]
type = SIMPLE
momentum_l_abs_tol = 1e-10
pressure_l_abs_tol = 1e-10
energy_l_abs_tol = 1e-10
momentum_l_tol = 0
pressure_l_tol = 0
energy_l_tol = 0
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
energy_system = 'energy_system'
momentum_equation_relaxation = 0.8
energy_equation_relaxation = 0.9
pressure_variable_relaxation = 0.3
num_iterations = 200
pressure_absolute_tolerance = 1e-10
momentum_absolute_tolerance = 1e-10
energy_absolute_tolerance = 1e-10
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
energy_petsc_options_iname = '-pc_type -pc_hypre_type'
energy_petsc_options_value = 'hypre boomeramg'
print_fields = false
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/steady-transient-compare/common-blocks.i)
mu = 2.6
rho = 1.0
advected_interp_method = 'average'
cp = 300
k = 20
alpha_b = 1e-4
[Mesh]
[mesh]
type = CartesianMeshGenerator
dim = 2
dx = '1.'
dy = '0.2'
ix = '10'
iy = '5'
[]
[]
[Problem]
linear_sys_names = 'u_system v_system pressure_system energy_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = RhieChowMassFlux
u = vel_x
v = vel_y
pressure = pressure
rho = ${rho}
p_diffusion_kernel = p_diffusion
[]
[]
[Variables]
[vel_x]
type = MooseLinearVariableFVReal
initial_condition = 0.5
solver_sys = u_system
[]
[vel_y]
type = MooseLinearVariableFVReal
solver_sys = v_system
initial_condition = 0.0
[]
[pressure]
type = MooseLinearVariableFVReal
solver_sys = pressure_system
initial_condition = 0.2
[]
[T]
type = MooseLinearVariableFVReal
solver_sys = energy_system
initial_condition = 300
[]
[]
[LinearFVKernels]
[u_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_x
advected_interp_method = ${advected_interp_method}
mu = ${mu}
u = vel_x
v = vel_y
momentum_component = 'x'
rhie_chow_user_object = 'rc'
use_nonorthogonal_correction = false
[]
[v_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_y
advected_interp_method = ${advected_interp_method}
mu = ${mu}
u = vel_x
v = vel_y
momentum_component = 'y'
rhie_chow_user_object = 'rc'
use_nonorthogonal_correction = false
[]
[u_pressure]
type = LinearFVMomentumPressure
variable = vel_x
pressure = pressure
momentum_component = 'x'
[]
[v_pressure]
type = LinearFVMomentumPressure
variable = vel_y
pressure = pressure
momentum_component = 'y'
[]
[u_boussinesq]
type = LinearFVMomentumBoussinesq
variable = vel_x
rho = ${rho}
gravity = '0 -9.81 0'
alpha_name = ${alpha_b}
ref_temperature = 300.0
T_fluid = T
momentum_component = 'x'
[]
[v_boussinesq]
type = LinearFVMomentumBoussinesq
variable = vel_y
rho = ${rho}
gravity = '0 -9.81 0'
alpha_name = ${alpha_b}
ref_temperature = 300.0
T_fluid = T
momentum_component = 'y'
[]
[p_diffusion]
type = LinearFVAnisotropicDiffusion
variable = pressure
diffusion_tensor = Ainv
use_nonorthogonal_correction = false
[]
[HbyA_divergence]
type = LinearFVDivergence
variable = pressure
face_flux = HbyA
force_boundary_execution = true
[]
[h_advection]
type = LinearFVEnergyAdvection
variable = T
advected_quantity = temperature
cp = ${cp}
advected_interp_method = ${advected_interp_method}
rhie_chow_user_object = 'rc'
[]
[conduction]
type = LinearFVDiffusion
variable = T
diffusion_coeff = ${k}
use_nonorthogonal_correction = false
[]
[]
[FunctorMaterials]
[constant_functors]
type = GenericFunctorMaterial
prop_names = 'cp alpha_b'
prop_values = '${cp} ${alpha_b}'
[]
[]
[LinearFVBCs]
[inlet-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_x
functor = '1.1'
[]
[inlet-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_y
functor = '0.0'
[]
[walls-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom'
variable = vel_x
functor = 0.0
[]
[walls-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom'
variable = vel_y
functor = 0.0
[]
[outlet_p]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'right'
variable = pressure
functor = 1.4
[]
[outlet_u]
type = LinearFVAdvectionDiffusionOutflowBC
variable = vel_x
use_two_term_expansion = false
boundary = right
[]
[outlet_v]
type = LinearFVAdvectionDiffusionOutflowBC
variable = vel_y
use_two_term_expansion = false
boundary = right
[]
[inlet_top_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T
functor = ${fparse 300.0}
boundary = 'left top'
[]
[bottom_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T
functor = ${fparse 350.0}
boundary = bottom
[]
[outlet_T]
type = LinearFVAdvectionDiffusionOutflowBC
variable = T
use_two_term_expansion = false
boundary = right
[]
[]
[Outputs]
exodus = true
execute_on = FINAL
[]
(modules/navier_stokes/test/tests/finite_volume/ins/natural_convection/linear_segregated/2d/diff_heated_cavity_linear_segregated.i)
################################################################################
# MATERIAL PROPERTIES
################################################################################
rho = 3279.
T_0 = 875.0
mu = 1.
k_cond = 38.0
cp = 640.
alpha = 3.26e-4
walls = 'right left top bottom'
[GlobalParams]
rhie_chow_user_object = 'ins_rhie_chow_interpolator'
advected_interp_method = 'upwind'
u = superficial_vel_x
v = superficial_vel_y
[]
[Problem]
linear_sys_names = 'u_system v_system pressure_system energy_system'
previous_nl_solution_required = true
[]
################################################################################
# GEOMETRY
################################################################################
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 30
ny = 30
[]
[]
################################################################################
# EQUATIONS: VARIABLES, KERNELS & BCS
################################################################################
[UserObjects]
[ins_rhie_chow_interpolator]
type = RhieChowMassFlux
u = superficial_vel_x
v = superficial_vel_y
pressure = pressure
rho = ${rho}
p_diffusion_kernel = p_diffusion
[]
[]
[Variables]
[superficial_vel_x]
type = MooseLinearVariableFVReal
solver_sys = u_system
[]
[superficial_vel_y]
type = MooseLinearVariableFVReal
solver_sys = v_system
[]
[pressure]
type = MooseLinearVariableFVReal
initial_condition = 0
solver_sys = pressure_system
[]
[T_fluid]
type = MooseLinearVariableFVReal
solver_sys = energy_system
initial_condition = 875
[]
[]
[LinearFVKernels]
[u_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = superficial_vel_x
mu = ${mu}
momentum_component = 'x'
use_nonorthogonal_correction = true
[]
[u_pressure]
type = LinearFVMomentumPressure
variable = superficial_vel_x
pressure = pressure
momentum_component = 'x'
[]
[u_buoyancy]
type = LinearFVMomentumBoussinesq
variable = superficial_vel_x
T_fluid = T_fluid
gravity = '0 -9.81 0'
rho = ${rho}
ref_temperature = ${T_0}
alpha_name = ${alpha}
momentum_component = 'x'
[]
[v_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = superficial_vel_y
mu = ${mu}
momentum_component = 'y'
use_nonorthogonal_correction = true
[]
[v_pressure]
type = LinearFVMomentumPressure
variable = superficial_vel_y
pressure = pressure
momentum_component = 'y'
[]
[v_buoyancy]
type = LinearFVMomentumBoussinesq
variable = superficial_vel_y
T_fluid = T_fluid
gravity = '0 -9.81 0'
rho = ${rho}
ref_temperature = ${T_0}
alpha_name = ${alpha}
momentum_component = 'y'
[]
[p_diffusion]
type = LinearFVAnisotropicDiffusion
variable = pressure
diffusion_tensor = Ainv
use_nonorthogonal_correction = true
[]
[HbyA_divergence]
type = LinearFVDivergence
variable = pressure
face_flux = HbyA
force_boundary_execution = true
[]
####### FUEL ENERGY EQUATION #######
[heat_advection]
type = LinearFVEnergyAdvection
variable = T_fluid
advected_quantity = temperature
cp = ${cp}
[]
[conduction]
type = LinearFVDiffusion
variable = T_fluid
diffusion_coeff = ${fparse k_cond}
[]
[]
[LinearFVBCs]
[no-slip-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = superficial_vel_x
boundary = ${walls}
functor = 0
[]
[no-slip-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = superficial_vel_y
boundary = ${walls}
functor = 0
[]
[T_cold]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T_fluid
boundary = 'right'
functor = 870.0
[]
[T_hot]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T_fluid
boundary = 'left'
functor = 880.0
[]
[T_all]
type = LinearFVAdvectionDiffusionExtrapolatedBC
variable = T_fluid
boundary = 'top bottom'
use_two_term_expansion = false
[]
[pressure-extrapolation]
type = LinearFVExtrapolatedPressureBC
boundary = ${walls}
variable = pressure
use_two_term_expansion = false
[]
[]
[FunctorMaterials]
[constant_functors]
type = GenericFunctorMaterial
prop_names = 'cp alpha_b'
prop_values = '${cp} ${alpha}'
[]
[]
################################################################################
# EXECUTION / SOLVE
################################################################################
[Executioner]
type = SIMPLE
momentum_l_abs_tol = 1e-11
pressure_l_abs_tol = 1e-11
energy_l_abs_tol = 1e-11
momentum_l_tol = 0
pressure_l_tol = 0
energy_l_tol = 0
rhie_chow_user_object = 'ins_rhie_chow_interpolator'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
energy_system = 'energy_system'
momentum_equation_relaxation = 0.3
pressure_variable_relaxation = 0.7
energy_equation_relaxation = 0.95
num_iterations = 3000
pressure_absolute_tolerance = 1e-8
momentum_absolute_tolerance = 1e-8
energy_absolute_tolerance = 1e-8
print_fields = false
momentum_l_max_its = 300
pin_pressure = true
pressure_pin_value = 0.0
pressure_pin_point = '0.5 0.0 0.0'
# momentum_petsc_options = '-ksp_monitor'
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
energy_petsc_options_iname = '-pc_type -pc_hypre_type'
energy_petsc_options_value = 'hypre boomeramg'
[]
################################################################################
# SIMULATION OUTPUTS
################################################################################
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/2d/2d-boussinesq-transient.i)
mu = 2.6
rho = 1.0
advected_interp_method = 'average'
cp = 300
k = 10
alpha_b = 1e-4
[Mesh]
[mesh]
type = CartesianMeshGenerator
dim = 2
dx = '1.'
dy = '0.2'
ix = '10'
iy = '5'
[]
[]
[Problem]
linear_sys_names = 'u_system v_system pressure_system energy_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = RhieChowMassFlux
u = vel_x
v = vel_y
pressure = pressure
rho = ${rho}
p_diffusion_kernel = p_diffusion
[]
[]
[Variables]
[vel_x]
type = MooseLinearVariableFVReal
initial_condition = 0.5
solver_sys = u_system
[]
[vel_y]
type = MooseLinearVariableFVReal
solver_sys = v_system
initial_condition = 0.0
[]
[pressure]
type = MooseLinearVariableFVReal
solver_sys = pressure_system
initial_condition = 0.2
[]
[T]
type = MooseLinearVariableFVReal
solver_sys = energy_system
initial_condition = 300
[]
[]
[LinearFVKernels]
[u_time]
type = LinearFVTimeDerivative
variable = vel_x
factor = ${rho}
[]
[v_time]
type = LinearFVTimeDerivative
variable = vel_y
factor = ${rho}
[]
[u_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_x
advected_interp_method = ${advected_interp_method}
mu = ${mu}
u = vel_x
v = vel_y
momentum_component = 'x'
rhie_chow_user_object = 'rc'
use_nonorthogonal_correction = false
[]
[v_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_y
advected_interp_method = ${advected_interp_method}
mu = ${mu}
u = vel_x
v = vel_y
momentum_component = 'y'
rhie_chow_user_object = 'rc'
use_nonorthogonal_correction = false
[]
[u_pressure]
type = LinearFVMomentumPressure
variable = vel_x
pressure = pressure
momentum_component = 'x'
[]
[v_pressure]
type = LinearFVMomentumPressure
variable = vel_y
pressure = pressure
momentum_component = 'y'
[]
[u_boussinesq]
type = LinearFVMomentumBoussinesq
variable = vel_x
rho = ${rho}
gravity = '0 -9.81 0'
alpha_name = ${alpha_b}
ref_temperature = 300.0
T_fluid = T
momentum_component = 'x'
[]
[v_boussinesq]
type = LinearFVMomentumBoussinesq
variable = vel_y
rho = ${rho}
gravity = '0 -9.81 0'
alpha_name = ${alpha_b}
ref_temperature = 300.0
T_fluid = T
momentum_component = 'y'
[]
[p_diffusion]
type = LinearFVAnisotropicDiffusion
variable = pressure
diffusion_tensor = Ainv
use_nonorthogonal_correction = false
[]
[HbyA_divergence]
type = LinearFVDivergence
variable = pressure
face_flux = HbyA
force_boundary_execution = true
[]
[h_time]
type = LinearFVTimeDerivative
variable = T
factor = ${fparse rho*cp}
[]
[h_advection]
type = LinearFVEnergyAdvection
variable = T
advected_quantity = temperature
cp = ${cp}
advected_interp_method = ${advected_interp_method}
rhie_chow_user_object = 'rc'
[]
[conduction]
type = LinearFVDiffusion
variable = T
diffusion_coeff = ${k}
use_nonorthogonal_correction = false
[]
[]
[FunctorMaterials]
[constant_functors]
type = GenericFunctorMaterial
prop_names = 'cp alpha_b'
prop_values = '${cp} ${alpha_b}'
[]
[]
[LinearFVBCs]
[inlet-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_x
functor = '1.1'
[]
[inlet-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_y
functor = '0.0'
[]
[walls-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom'
variable = vel_x
functor = 0.0
[]
[walls-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom'
variable = vel_y
functor = 0.0
[]
[outlet_p]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'right'
variable = pressure
functor = 1.4
[]
[outlet_u]
type = LinearFVAdvectionDiffusionOutflowBC
variable = vel_x
use_two_term_expansion = false
boundary = right
[]
[outlet_v]
type = LinearFVAdvectionDiffusionOutflowBC
variable = vel_y
use_two_term_expansion = false
boundary = right
[]
[inlet_top_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T
functor = 300.0
boundary = 'left top'
[]
[bottom_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T
functor = wall-temperature
boundary = bottom
[]
[outlet_T]
type = LinearFVAdvectionDiffusionOutflowBC
variable = T
use_two_term_expansion = false
boundary = right
[]
[]
[Functions]
[wall-temperature]
type = ParsedFunction
expression = '350 + 50 * sin(6.28*t)'
[]
[]
[Executioner]
type = PIMPLE
momentum_l_abs_tol = 1e-12
pressure_l_abs_tol = 1e-12
energy_l_abs_tol = 1e-12
momentum_l_tol = 1e-12
pressure_l_tol = 1e-12
energy_l_tol = 1e-12
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
energy_system = 'energy_system'
momentum_equation_relaxation = 0.8
pressure_variable_relaxation = 0.3
energy_equation_relaxation = 0.9
num_iterations = 100
pressure_absolute_tolerance = 1e-11
momentum_absolute_tolerance = 1e-11
energy_absolute_tolerance = 1e-11
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
energy_petsc_options_iname = '-pc_type -pc_hypre_type'
energy_petsc_options_value = 'hypre boomeramg'
print_fields = false
continue_on_max_its = true
dt = 0.01
num_steps = 6
num_piso_iterations = 0
[]
[Outputs]
exodus = true
[]