- advective_flux_calculatorAdvectiveFluxCalculator UserObject
C++ Type:UserObjectName
Description:AdvectiveFluxCalculator UserObject
- variableThe name of the variable that this Kernel operates on
C++ Type:NonlinearVariableName
Description:The name of the variable that this Kernel operates on
FluxLimitedTVDAdvection
This Kernel
implements the weak form of It employs Kuzmin-Turek (Kuzmin and Turek, 2004) stabilization, so needs a corresponding AdvectiveFluxCalculatorConstantVelocity UserObject. This sort of stabilization is described in detail in A worked example of Kuzmin-Turek stabilization. Also see numerical diffusion for example of how the Kuzmin-Turek scheme compares with other numerical schemes.
Input Parameters
- blockThe list of block ids (SubdomainID) that this object will be applied
C++ Type:std::vector
Options:
Description:The list of block ids (SubdomainID) that this object will be applied
- displacementsThe displacements
C++ Type:std::vector
Options:
Description:The displacements
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector
Options:
Description:Adds user-defined labels for accessing object parameters via control logic.
- diag_save_inThe name of auxiliary variables to save this Kernel's diagonal Jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector
Options:
Description:The name of auxiliary variables to save this Kernel's diagonal Jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Options:
Description:Set the enabled status of the MooseObject.
- implicitTrueDetermines whether this object is calculated using an implicit or explicit form
Default:True
C++ Type:bool
Options:
Description:Determines whether this object is calculated using an implicit or explicit form
- save_inThe name of auxiliary variables to save this Kernel's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector
Options:
Description:The name of auxiliary variables to save this Kernel's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
- seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Options:
Description:The seed for the master random number generator
- 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
- extra_matrix_tagsThe extra tags for the matrices this Kernel should fill
C++ Type:std::vector
Options:
Description:The extra tags for the matrices this Kernel should fill
- extra_vector_tagsThe extra tags for the vectors this Kernel should fill
C++ Type:std::vector
Options:
Description:The extra tags for the vectors this Kernel should fill
- matrix_tagssystemThe tag for the matrices this Kernel should fill
Default:system
C++ Type:MultiMooseEnum
Options:nontime system
Description:The tag for the matrices this Kernel should fill
- vector_tagsnontimeThe tag for the vectors this Kernel should fill
Default:nontime
C++ Type:MultiMooseEnum
Options:nontime time
Description:The tag for the vectors this Kernel should fill
Tagging Parameters
Input Files
- modules/porous_flow/test/tests/flux_limited_TVD_advection/jacobian_02.i
- modules/porous_flow/test/tests/flux_limited_TVD_advection/fltvd_1D_adaptivity.i
- modules/porous_flow/test/tests/flux_limited_TVD_advection/except_01.i
- modules/porous_flow/test/tests/flux_limited_TVD_advection/fltvd_2D_trimesh.i
- modules/porous_flow/test/tests/flux_limited_TVD_advection/fltvd_2D_angle.i
- modules/porous_flow/test/tests/flux_limited_TVD_advection/jacobian_01.i
- modules/porous_flow/test/tests/numerical_diffusion/fltvd_no_antidiffusion.i
- modules/porous_flow/test/tests/flux_limited_TVD_advection/jacobian_03.i
- modules/porous_flow/test/tests/flux_limited_TVD_advection/fltvd_2D_blocks.i
- modules/porous_flow/test/tests/flux_limited_TVD_advection/fltvd_2D.i
- modules/porous_flow/test/tests/flux_limited_TVD_advection/fltvd_1D.i
- modules/porous_flow/test/tests/flux_limited_TVD_advection/fltvd_3D.i
- modules/porous_flow/test/tests/numerical_diffusion/fltvd.i
- modules/porous_flow/test/tests/numerical_diffusion/fltvd_none.i
modules/porous_flow/test/tests/flux_limited_TVD_advection/jacobian_02.i
# Checking the Jacobian of Flux-Limited TVD Advection, using flux_limiter_type = superbee
# Here we use snes_check_jacobian instead of snes_type=test. The former just checks the Jacobian for the
# random initial conditions, while the latter checks for u=1 and u=-1
#
# The Jacobian is correct for u=1 and u=-1, but the finite-difference scheme used by snes_type=test gives the
# wrong answer.
# For u=1, the Kuzmin-Turek scheme adds as much antidiffusion as possible, resulting in a central-difference
# version of advection (flux_limiter = 1). This is correct, and the Jacobian is calculated correctly.
# However, when computing the Jacobian using finite differences, u is increased or decreased at a node.
# This results in that node being at a maximum or minimum, which means no antidiffusion should be added
# (flux_limiter = 0). This corresponds to a full-upwind scheme. So the finite-difference computes the
# Jacobian in the full-upwind scenario, which is incorrect (the original residual = 0, after finite-differencing
# the residual comes from the full-upwind scenario).
[Mesh]
type = GeneratedMesh
dim = 3
nx = 2
xmin = 0
xmax = 1
ny = 2
ymin = -1
ymax = 2
bias_y = 1.5
nz = 2
zmin = 1
zmax = 2
bias_z = 0.8
[]
[Variables]
[./u]
[../]
[]
[ICs]
[./u]
type = FunctionIC
variable = u
function = 'x + 0.5 * y - 0.4 * z - 0.1 * sin(x) - 0.1 * cos(y) + 0.2 * exp(-z)'
[../]
[]
[Kernels]
[./flux]
type = FluxLimitedTVDAdvection
variable = u
advective_flux_calculator = fluo
[../]
[]
[UserObjects]
[./fluo]
type = AdvectiveFluxCalculatorConstantVelocity
flux_limiter_type = superbee
u = u
velocity = '1 -2 1.5'
[../]
[]
[Preconditioning]
active = smp
[./smp]
type = SMP
full = true
petsc_options = '-snes_check_jacobian'
[../]
[]
[Executioner]
type = Transient
solve_type = Linear # this is to force convergence even though the nonlinear residual is high: we just care about the Jacobian in this test
end_time = 1
num_steps = 1
dt = 1
[]
modules/porous_flow/test/tests/flux_limited_TVD_advection/fltvd_1D_adaptivity.i
# Using Flux-Limited TVD Advection ala Kuzmin and Turek, with antidiffusion from superbee flux limiting
# 1D version with mesh adaptivity
[Mesh]
type = GeneratedMesh
dim = 1
nx = 10
xmin = 0
xmax = 1
[]
[Adaptivity]
initial_steps = 1
initial_marker = tracer_marker
marker = tracer_marker
max_h_level = 1
[./Markers]
[./tracer_marker]
type = ValueRangeMarker
variable = tracer
lower_bound = 0.02
upper_bound = 0.98
[../]
[../]
[]
[Variables]
[./tracer]
[../]
[]
[ICs]
[./tracer]
type = FunctionIC
variable = tracer
function = 'if(x<0.1,0,if(x>0.3,0,1))'
[../]
[]
[Kernels]
[./mass_dot]
type = MassLumpedTimeDerivative
variable = tracer
[../]
[./flux]
type = FluxLimitedTVDAdvection
variable = tracer
advective_flux_calculator = fluo
[../]
[]
[UserObjects]
[./fluo]
type = AdvectiveFluxCalculatorConstantVelocity
flux_limiter_type = superbee
u = tracer
velocity = '0.1 0 0'
[../]
[]
[BCs]
[./no_tracer_on_left]
type = DirichletBC
variable = tracer
value = 0
boundary = left
[../]
[./remove_tracer]
# Ideally, an OutflowBC would be used, but that does not exist in the framework
# In 1D VacuumBC is the same as OutflowBC, with the alpha parameter being twice the velocity
type = VacuumBC
boundary = right
alpha = 0.2 # 2 * velocity
variable = tracer
[../]
[]
[Preconditioning]
active = basic
[./basic]
type = SMP
full = true
petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
petsc_options_value = ' asm lu NONZERO 2'
[../]
[./preferred_but_might_not_be_installed]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = ' lu mumps'
[../]
[]
[VectorPostprocessors]
[./tracer]
type = LineValueSampler
start_point = '0 0 0'
end_point = '1 0 0'
num_points = 11
sort_by = x
variable = tracer
[../]
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 6
dt = 6E-2
nl_abs_tol = 1E-8
nl_max_its = 500
timestep_tolerance = 1E-3
[]
[Outputs]
print_linear_residuals = false
[./out]
type = CSV
execute_on = final
[../]
[]
modules/porous_flow/test/tests/flux_limited_TVD_advection/except_01.i
# Exception test that AdvectiveFluxCalculator is indeed executed on linear
[Mesh]
type = GeneratedMesh
dim = 1
nx = 1
[]
[Variables]
[./u]
[../]
[]
[Kernels]
[./flux]
type = FluxLimitedTVDAdvection
variable = u
advective_flux_calculator = fluo
[../]
[]
[UserObjects]
[./fluo]
type = AdvectiveFluxCalculatorConstantVelocity
execute_on = 'nonlinear timestep_begin timestep_end final initial'
u = u
velocity = '0 0 0'
[../]
[]
[Preconditioning]
active = smp
[./smp]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 1
num_steps = 1
dt = 1
[]
modules/porous_flow/test/tests/flux_limited_TVD_advection/fltvd_2D_trimesh.i
# Using Flux-Limited TVD Advection ala Kuzmin and Turek, with antidiffusion from superbee flux limiting
# 2D version
[Mesh]
type = FileMesh
file = trimesh.msh
[]
[GlobalParams]
block = '50'
[]
[Variables]
[./tracer]
[../]
[]
[ICs]
[./tracer]
type = FunctionIC
variable = tracer
function = 'if(x<0.1,0,if(x>0.305,0,1))'
[../]
[]
[Kernels]
[./mass_dot]
type = MassLumpedTimeDerivative
variable = tracer
[../]
[./flux]
type = FluxLimitedTVDAdvection
variable = tracer
advective_flux_calculator = fluo
[../]
[]
[UserObjects]
[./fluo]
type = AdvectiveFluxCalculatorConstantVelocity
flux_limiter_type = superbee
u = tracer
velocity = '0.1 0 0'
[../]
[]
[BCs]
[./no_tracer_on_left]
type = DirichletBC
variable = tracer
value = 0
boundary = left
[../]
[./remove_tracer]
# Ideally, an OutflowBC would be used, but that does not exist in the framework
# In 1D VacuumBC is the same as OutflowBC, with the alpha parameter being twice the velocity
type = VacuumBC
boundary = right
alpha = 0.2 # 2 * velocity
variable = tracer
[../]
[]
[Preconditioning]
active = basic
[./basic]
type = SMP
full = true
petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
petsc_options_value = ' asm lu NONZERO 2'
[../]
[./preferred_but_might_not_be_installed]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = ' lu mumps'
[../]
[]
[VectorPostprocessors]
[./tracer]
type = LineValueSampler
start_point = '0 0 0'
end_point = '1 0.04 0'
num_points = 101
sort_by = x
variable = tracer
[../]
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 6
dt = 6E-2
timestep_tolerance = 1E-3
[]
[Outputs]
print_linear_residuals = false
[./out]
type = CSV
execute_on = final
[../]
[]
modules/porous_flow/test/tests/flux_limited_TVD_advection/fltvd_2D_angle.i
# Using Flux-Limited TVD Advection ala Kuzmin and Turek, with antidiffusion from superbee flux limiting
# 2D version with velocity = (0.1, 0.2, 0)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
xmin = 0
xmax = 1
ny = 10
ymin = 0
ymax = 1
[]
[Variables]
[./tracer]
[../]
[]
[ICs]
[./tracer]
type = FunctionIC
variable = tracer
function = 'if(x<0.1 | x > 0.3 | y < 0.1 | y > 0.3, 0, 1)'
[../]
[]
[Kernels]
[./mass_dot]
type = MassLumpedTimeDerivative
variable = tracer
[../]
[./flux]
type = FluxLimitedTVDAdvection
variable = tracer
advective_flux_calculator = fluo
[../]
[]
[UserObjects]
[./fluo]
type = AdvectiveFluxCalculatorConstantVelocity
flux_limiter_type = superbee
u = tracer
velocity = '0.1 0.2 0'
[../]
[]
[BCs]
[./no_tracer_on_left]
type = DirichletBC
variable = tracer
value = 0
boundary = left
[../]
[./remove_tracer]
# Ideally, an OutflowBC would be used, but that does not exist in the framework
# In 1D VacuumBC is the same as OutflowBC, with the alpha parameter being twice the velocity
type = VacuumBC
boundary = right
alpha = 0.2 # 2 * velocity
variable = tracer
[../]
[]
[Preconditioning]
active = basic
[./basic]
type = SMP
full = true
petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
petsc_options_value = ' asm lu NONZERO 2'
[../]
[./preferred_but_might_not_be_installed]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = ' lu mumps'
[../]
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 2
dt = 0.1
[]
[Outputs]
print_linear_residuals = false
[./out]
type = Exodus
execute_on = 'initial final'
[../]
[]
modules/porous_flow/test/tests/flux_limited_TVD_advection/jacobian_01.i
# Checking the Jacobian of Flux-Limited TVD Advection, using flux_limiter_type = none
[Mesh]
type = GeneratedMesh
dim = 3
nx = 3
xmin = 0
xmax = 1
ny = 4
ymin = -1
ymax = 2
bias_y = 1.5
nz = 4
zmin = 1
zmax = 2
bias_z = 0.8
[]
[Variables]
[./u]
[../]
[]
[ICs]
[./u]
type = RandomIC
variable = u
[../]
[]
[Kernels]
[./flux]
type = FluxLimitedTVDAdvection
variable = u
advective_flux_calculator = fluo
[../]
[]
[UserObjects]
[./fluo]
type = AdvectiveFluxCalculatorConstantVelocity
flux_limiter_type = none
u = u
velocity = '1 -2 1.5'
[../]
[]
[Preconditioning]
active = smp
[./smp]
type = SMP
full = true
petsc_options_iname = '-snes_type'
petsc_options_value = 'test'
[../]
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 1
num_steps = 1
dt = 1
[]
modules/porous_flow/test/tests/numerical_diffusion/fltvd_no_antidiffusion.i
# Using Flux-Limited TVD Advection ala Kuzmin and Turek, but without any antidiffusion
[Mesh]
type = GeneratedMesh
dim = 1
nx = 100
xmin = 0
xmax = 1
[]
[Variables]
[./tracer]
[../]
[]
[ICs]
[./tracer]
type = FunctionIC
variable = tracer
function = 'if(x<0.1,0,if(x>0.3,0,1))'
[../]
[]
[Kernels]
[./mass_dot]
type = MassLumpedTimeDerivative
variable = tracer
[../]
[./flux]
type = FluxLimitedTVDAdvection
variable = tracer
advective_flux_calculator = fluo
[../]
[]
[UserObjects]
[./fluo]
type = AdvectiveFluxCalculatorConstantVelocity
flux_limiter_type = none
u = tracer
velocity = '0.1 0 0'
[../]
[]
[BCs]
[./no_tracer_on_left]
type = DirichletBC
variable = tracer
value = 0
boundary = left
[../]
[./remove_tracer]
# Ideally, an OutflowBC would be used, but that does not exist in the framework
# In 1D VacuumBC is the same as OutflowBC, with the alpha parameter being twice the velocity
type = VacuumBC
boundary = right
alpha = 0.2 # 2 * velocity
variable = tracer
[../]
[]
[Preconditioning]
active = basic
[./basic]
type = SMP
full = true
petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
petsc_options_value = ' asm lu NONZERO 2'
[../]
[./preferred_but_might_not_be_installed]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = ' lu mumps'
[../]
[]
[VectorPostprocessors]
[./tracer]
type = LineValueSampler
start_point = '0 0 0'
end_point = '1 0 0'
num_points = 101
sort_by = x
variable = tracer
[../]
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 6
dt = 6E-1
nl_abs_tol = 1E-8
nl_max_its = 500
timestep_tolerance = 1E-3
[]
[Outputs]
[./out]
type = CSV
execute_on = final
[../]
[]
modules/porous_flow/test/tests/flux_limited_TVD_advection/jacobian_03.i
# Checking the Jacobian of Flux-Limited TVD Advection, using flux_limiter_type = vanleer
#
# The initial conditions are u=x. This means that the argument of the flux limiter is 1, so that
# the flux_limiter=1 everywhere, irrespective of flux_limiter_type (except for 'none'). However
# superbee and minmod are nondifferentiable at this point, so using those flux_limiter_type will
# result in a poor Jacobian
#
# Here we use snes_check_jacobian instead of snes_type=test. The former just checks the Jacobian for the
# random initial conditions, while the latter checks for u=1 and u=-1
#
# The Jacobian is correct for u=1 and u=-1, but the finite-difference scheme used by snes_type=test gives the
# wrong answer.
# For u=1, the Kuzmin-Turek scheme adds as much antidiffusion as possible, resulting in a central-difference
# version of advection (flux_limiter = 1). This is correct, and the Jacobian is calculated correctly.
# However, when computing the Jacobian using finite differences, u is increased or decreased at a node.
# This results in that node being at a maximum or minimum, which means no antidiffusion should be added
# (flux_limiter = 0). This corresponds to a full-upwind scheme. So the finite-difference computes the
# Jacobian in the full-upwind scenario, which is incorrect (the original residual = 0, after finite-differencing
# the residual comes from the full-upwind scenario).
[Mesh]
type = GeneratedMesh
dim = 1
nx = 6
[]
[Variables]
[./u]
[../]
[]
[ICs]
[./u]
type = FunctionIC
variable = u
function = 'x'
[../]
[]
[Kernels]
[./flux]
type = FluxLimitedTVDAdvection
variable = u
advective_flux_calculator = fluo
[../]
[]
[UserObjects]
[./fluo]
type = AdvectiveFluxCalculatorConstantVelocity
flux_limiter_type = vanleer
u = u
velocity = '1 -2 1.5'
[../]
[]
[Preconditioning]
active = smp
[./smp]
type = SMP
full = true
petsc_options = '-snes_check_jacobian'
[../]
[]
[Executioner]
type = Transient
solve_type = Linear # this is to force convergence even though the nonlinear residual is high: we just care about the Jacobian in this test
end_time = 1
num_steps = 1
dt = 1
[]
modules/porous_flow/test/tests/flux_limited_TVD_advection/fltvd_2D_blocks.i
# Using Flux-Limited TVD Advection ala Kuzmin and Turek
# 2D version with blocks
# Top block: tracer is defined here, with velocity = (0.1, 0, 0)
# Central block: tracer is not defined here
# Bottom block: tracer is defined here, with velocity = (-0.1, 0, 0)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 10
xmin = 0
xmax = 1
ny = 5
ymin = 0
ymax = 1
[]
[./top]
input = gen
type = SubdomainBoundingBoxGenerator
bottom_left = '0 0.6 0'
top_right = '1 1 0'
block_id = 1
[../]
[./center]
input = bottom
type = SubdomainBoundingBoxGenerator
bottom_left = '0 0.4 0'
top_right = '1 0.6 0'
block_id = 2
[../]
[./bottom]
input = top
type = SubdomainBoundingBoxGenerator
bottom_left = '0 0 0'
top_right = '1 0.6 0'
block_id = 3
[../]
[./split_bdys]
type = BreakBoundaryOnSubdomainGenerator
input = center
boundaries = 'left right'
[../]
[]
[GlobalParams]
block = '1 2 3'
[]
[Variables]
[./tracer]
block = '1 3'
[../]
[./dummy]
[../]
[]
[ICs]
[./tracer_top]
type = FunctionIC
variable = tracer
function = 'if(x<0.1 | x>0.3, 0, 1)'
block = '1'
[../]
[./tracer_bot]
type = FunctionIC
variable = tracer
function = 'if(x<0.7 | x > 0.9, 0, 1)'
block = '3'
[../]
[]
[Kernels]
[./mass_dot]
type = MassLumpedTimeDerivative
variable = tracer
block = '1 3'
[../]
[./flux_top]
type = FluxLimitedTVDAdvection
variable = tracer
advective_flux_calculator = fluo_top
block = '1'
[../]
[./flux_bot]
type = FluxLimitedTVDAdvection
variable = tracer
advective_flux_calculator = fluo_bot
block = '3'
[../]
[.dummy]
type = TimeDerivative
variable = dummy
[../]
[]
[UserObjects]
[./fluo_top]
type = AdvectiveFluxCalculatorConstantVelocity
flux_limiter_type = superbee
u = tracer
velocity = '0.1 0 0'
block = '1'
[../]
[./fluo_bot]
type = AdvectiveFluxCalculatorConstantVelocity
flux_limiter_type = superbee
u = tracer
velocity = '-0.1 0 0'
block = '3'
[../]
[]
[BCs]
[./no_tracer_on_left_top]
type = DirichletBC
variable = tracer
value = 0
boundary = 'left_to_1'
[../]
[./remove_tracer_top]
# Ideally, an OutflowBC would be used, but that does not exist in the framework
# In 1D VacuumBC is the same as OutflowBC, with the alpha parameter being twice the velocity
type = VacuumBC
boundary = 'right_to_1'
alpha = 0.2 # 2 * velocity
variable = tracer
[../]
[./no_tracer_on_left_bot]
# Ideally, an OutflowBC would be used, but that does not exist in the framework
# In 1D VacuumBC is the same as OutflowBC, with the alpha parameter being twice the velocity
type = VacuumBC
boundary = 'left_to_3'
alpha = 0.2 # 2 * velocity
variable = tracer
[../]
[./remove_tracer_bot]
type = DirichletBC
variable = tracer
value = 0
boundary = 'right_to_3'
[../]
[]
[Preconditioning]
active = basic
[./basic]
type = SMP
full = true
petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
petsc_options_value = ' asm lu NONZERO 2'
[../]
[./preferred_but_might_not_be_installed]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = ' lu mumps'
[../]
[]
[VectorPostprocessors]
[./tracer_bot]
type = LineValueSampler
start_point = '0 0 0'
end_point = '1 0 0'
num_points = 11
sort_by = x
variable = tracer
[../]
[./tracer_top]
type = LineValueSampler
start_point = '0 1 0'
end_point = '1 1 0'
num_points = 11
sort_by = x
variable = tracer
[../]
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 6
dt = 6E-2
timestep_tolerance = 1E-3
[]
[Outputs]
print_linear_residuals = false
[./out]
type = CSV
execute_on = final
[../]
[]
modules/porous_flow/test/tests/flux_limited_TVD_advection/fltvd_2D.i
# Using Flux-Limited TVD Advection ala Kuzmin and Turek, with antidiffusion from superbee flux limiting
# 2D version
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
xmin = 0
xmax = 1
ny = 4
ymin = 0
ymax = 0.5
[]
[Variables]
[./tracer]
[../]
[]
[ICs]
[./tracer]
type = FunctionIC
variable = tracer
function = 'if(x<0.1,0,if(x>0.3,0,1))'
[../]
[]
[Kernels]
[./mass_dot]
type = MassLumpedTimeDerivative
variable = tracer
[../]
[./flux]
type = FluxLimitedTVDAdvection
variable = tracer
advective_flux_calculator = fluo
[../]
[]
[UserObjects]
[./fluo]
type = AdvectiveFluxCalculatorConstantVelocity
flux_limiter_type = superbee
u = tracer
velocity = '0.1 0 0'
[../]
[]
[BCs]
[./no_tracer_on_left]
type = DirichletBC
variable = tracer
value = 0
boundary = left
[../]
[./remove_tracer]
# Ideally, an OutflowBC would be used, but that does not exist in the framework
# In 1D VacuumBC is the same as OutflowBC, with the alpha parameter being twice the velocity
type = VacuumBC
boundary = right
alpha = 0.2 # 2 * velocity
variable = tracer
[../]
[]
[Preconditioning]
active = basic
[./basic]
type = SMP
full = true
petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
petsc_options_value = ' asm lu NONZERO 2'
[../]
[./preferred_but_might_not_be_installed]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = ' lu mumps'
[../]
[]
[VectorPostprocessors]
[./tracer]
type = LineValueSampler
start_point = '0 0 0'
end_point = '1 0.5 0'
num_points = 11
sort_by = x
variable = tracer
[../]
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 6
dt = 6E-2
nl_abs_tol = 1E-8
nl_max_its = 500
timestep_tolerance = 1E-3
[]
[Outputs]
print_linear_residuals = false
[./out]
type = CSV
execute_on = final
[../]
[]
modules/porous_flow/test/tests/flux_limited_TVD_advection/fltvd_1D.i
# Using Flux-Limited TVD Advection ala Kuzmin and Turek, with antidiffusion from superbee flux limiting
# 1D version
[Mesh]
type = GeneratedMesh
dim = 1
nx = 10
xmin = 0
xmax = 1
[]
[Variables]
[./tracer]
[../]
[]
[ICs]
[./tracer]
type = FunctionIC
variable = tracer
function = 'if(x<0.1,0,if(x>0.3,0,1))'
[../]
[]
[Kernels]
[./mass_dot]
type = MassLumpedTimeDerivative
variable = tracer
[../]
[./flux]
type = FluxLimitedTVDAdvection
variable = tracer
advective_flux_calculator = fluo
[../]
[]
[UserObjects]
[./fluo]
type = AdvectiveFluxCalculatorConstantVelocity
flux_limiter_type = superbee
u = tracer
velocity = '0.1 0 0'
[../]
[]
[BCs]
[./no_tracer_on_left]
type = DirichletBC
variable = tracer
value = 0
boundary = left
[../]
[./remove_tracer]
# Ideally, an OutflowBC would be used, but that does not exist in the framework
# In 1D VacuumBC is the same as OutflowBC, with the alpha parameter being twice the velocity
type = VacuumBC
boundary = right
alpha = 0.2 # 2 * velocity
variable = tracer
[../]
[]
[Preconditioning]
active = basic
[./basic]
type = SMP
full = true
petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
petsc_options_value = ' asm lu NONZERO 2'
[../]
[./preferred_but_might_not_be_installed]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = ' lu mumps'
[../]
[]
[VectorPostprocessors]
[./tracer]
type = LineValueSampler
start_point = '0 0 0'
end_point = '1 0 0'
num_points = 11
sort_by = x
variable = tracer
[../]
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 6
dt = 6E-2
nl_abs_tol = 1E-8
nl_max_its = 500
timestep_tolerance = 1E-3
[]
[Outputs]
print_linear_residuals = false
[./out]
type = CSV
execute_on = final
[../]
[]
modules/porous_flow/test/tests/flux_limited_TVD_advection/fltvd_3D.i
# Using Flux-Limited TVD Advection ala Kuzmin and Turek, with antidiffusion from superbee flux limiting
# 3D version
[Mesh]
type = GeneratedMesh
dim = 3
nx = 10
xmin = 0
xmax = 1
ny = 4
ymin = 0
ymax = 0.5
nz = 3
zmin = 0
zmax = 2
[]
[Variables]
[./tracer]
[../]
[]
[Problem]
error_on_jacobian_nonzero_reallocation=true
[]
[ICs]
[./tracer]
type = FunctionIC
variable = tracer
function = 'if(x<0.1,0,if(x>0.3,0,1))'
[../]
[]
[Kernels]
[./mass_dot]
type = MassLumpedTimeDerivative
variable = tracer
[../]
[./flux]
type = FluxLimitedTVDAdvection
variable = tracer
advective_flux_calculator = fluo
[../]
[]
[UserObjects]
[./fluo]
type = AdvectiveFluxCalculatorConstantVelocity
flux_limiter_type = superbee
u = tracer
velocity = '0.1 0 0'
[../]
[]
[BCs]
[./no_tracer_on_left]
type = DirichletBC
variable = tracer
value = 0
boundary = left
[../]
[./remove_tracer]
# Ideally, an OutflowBC would be used, but that does not exist in the framework
# In 1D VacuumBC is the same as OutflowBC, with the alpha parameter being twice the velocity
type = VacuumBC
boundary = right
alpha = 0.2 # 2 * velocity
variable = tracer
[../]
[]
[Preconditioning]
active = basic
[./basic]
type = SMP
full = true
petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
petsc_options_value = ' asm lu NONZERO 2'
[../]
[./preferred_but_might_not_be_installed]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = ' lu mumps'
[../]
[]
[VectorPostprocessors]
[./tracer]
type = LineValueSampler
start_point = '0 0 0'
end_point = '1 0.5 2'
num_points = 11
sort_by = x
variable = tracer
[../]
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 6
dt = 6E-2
nl_abs_tol = 1E-8
nl_max_its = 500
timestep_tolerance = 1E-3
[]
[Outputs]
print_linear_residuals = false
[./out]
type = CSV
execute_on = final
[../]
[]
modules/porous_flow/test/tests/numerical_diffusion/fltvd.i
# Using Flux-Limited TVD Advection ala Kuzmin and Turek, with antidiffusion from superbee flux limiting
[Mesh]
type = GeneratedMesh
dim = 1
nx = 100
xmin = 0
xmax = 1
[]
[Variables]
[./tracer]
[../]
[]
[ICs]
[./tracer]
type = FunctionIC
variable = tracer
function = 'if(x<0.1,0,if(x>0.3,0,1))'
[../]
[]
[Kernels]
[./mass_dot]
type = MassLumpedTimeDerivative
variable = tracer
[../]
[./flux]
type = FluxLimitedTVDAdvection
variable = tracer
advective_flux_calculator = fluo
[../]
[]
[UserObjects]
[./fluo]
type = AdvectiveFluxCalculatorConstantVelocity
flux_limiter_type = superbee
u = tracer
velocity = '0.1 0 0'
[../]
[]
[BCs]
[./no_tracer_on_left]
type = DirichletBC
variable = tracer
value = 0
boundary = left
[../]
[./remove_tracer]
# Ideally, an OutflowBC would be used, but that does not exist in the framework
# In 1D VacuumBC is the same as OutflowBC, with the alpha parameter being twice the velocity
type = VacuumBC
boundary = right
alpha = 0.2 # 2 * velocity
variable = tracer
[../]
[]
[Preconditioning]
active = basic
[./basic]
type = SMP
full = true
petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
petsc_options_value = ' asm lu NONZERO 2'
[../]
[./preferred_but_might_not_be_installed]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = ' lu mumps'
[../]
[]
[VectorPostprocessors]
[./tracer]
type = LineValueSampler
start_point = '0 0 0'
end_point = '1 0 0'
num_points = 101
sort_by = x
variable = tracer
[../]
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 6
dt = 6E-2
nl_abs_tol = 1E-8
nl_max_its = 500
timestep_tolerance = 1E-3
[]
[Outputs]
[./out]
type = CSV
execute_on = final
[../]
[]
modules/porous_flow/test/tests/numerical_diffusion/fltvd_none.i
# Using Flux-Limited TVD Advection ala Kuzmin and Turek
# No antidiffusion, so this is identical to full-upwinding
[Mesh]
type = GeneratedMesh
dim = 1
nx = 100
xmin = 0
xmax = 1
[]
[Variables]
[./tracer]
[../]
[]
[ICs]
[./tracer]
type = FunctionIC
variable = tracer
function = 'if(x<0.1,0,if(x>0.3,0,1))'
[../]
[]
[Kernels]
[./mass_dot]
type = MassLumpedTimeDerivative
variable = tracer
[../]
[./flux]
type = FluxLimitedTVDAdvection
variable = tracer
advective_flux_calculator = fluo
[../]
[]
[UserObjects]
[./fluo]
type = AdvectiveFluxCalculator
flux_limiter_type = none
u = tracer
velocity = '0.1 0 0'
[../]
[]
[BCs]
[./no_tracer_on_left]
type = DirichletBC
variable = tracer
value = 0
boundary = left
[../]
[./remove_tracer]
# Ideally, an OutflowBC would be used, but that does not exist in the framework
# In 1D VacuumBC is the same as OutflowBC, with the alpha parameter being twice the velocity
type = VacuumBC
boundary = right
alpha = 0.2 # 2 * velocity
variable = tracer
[../]
[]
[Preconditioning]
active = basic
[./basic]
type = SMP
full = true
petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
petsc_options_value = ' asm lu NONZERO 2'
[../]
[./preferred_but_might_not_be_installed]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = ' lu mumps'
[../]
[]
[VectorPostprocessors]
[./tracer]
type = LineValueSampler
start_point = '0 0 0'
end_point = '1 0 0'
num_points = 101
sort_by = x
variable = tracer
[../]
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 6
dt = 6E-1
nl_abs_tol = 1E-8
nl_max_its = 500
timestep_tolerance = 1E-3
[]
[Outputs]
csv = true
execute_on = final
[]
- D. Kuzmin and S. Turek.
High-resolution FEM-TVD shcemes based on a fully multidimensional flux limiter.
Journal of Computational Physics, 198:131–158, 2004.[BibTeX]
@article{KuzminTurek2004, author = "Kuzmin, D. and Turek, S.", title = "{High-resolution FEM-TVD shcemes based on a fully multidimensional flux limiter}", journal = "Journal of Computational Physics", volume = "198", pages = "131--158", year = "2004" }