- arrayFalseTrue to make this variable a array variable regardless of number of components. If 'components' > 1, this will automatically be set to true.
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:True to make this variable a array variable regardless of number of components. If 'components' > 1, this will automatically be set to true.
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
Unit:(no unit assumed)
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
- components1Number of components for an array variable
Default:1
C++ Type:unsigned int
Unit:(no unit assumed)
Controllable:No
Description:Number of components for an array variable
- familyMONOMIALSpecifies the family of FE shape functions to use for this variable.
Default:MONOMIAL
C++ Type:MooseEnum
Unit:(no unit assumed)
Controllable:No
Description:Specifies the family of FE shape functions to use for this variable.
- fvTrueTrue to make this variable a finite volume variable
Default:True
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:True to make this variable a finite volume variable
- orderCONSTANTOrder of the FE shape function to use for this variable (additional orders not listed here are allowed, depending on the family).
Default:CONSTANT
C++ Type:MooseEnum
Unit:(no unit assumed)
Controllable:No
Description:Order of the FE shape function to use for this variable (additional orders not listed here are allowed, depending on the family).
- solver_sysnl0If this variable is a solver variable, this is the solver system to which it should be added.
Default:nl0
C++ Type:SolverSystemName
Unit:(no unit assumed)
Controllable:No
Description:If this variable is a solver variable, this is the solver system to which it should be added.
- use_dualFalseTrue to use dual basis for Lagrange multipliers
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:True to use dual basis for Lagrange multipliers
MooseVariableLinearFVReal
Overview
This variable is used for finite volume simulations where a linear system is solved with the following form:
where is the vector containing the degrees of freedom for the declared variable, while and are the corresponding system matrix and right hand side, respectively.
This variable can only be used with linear systems and cannot be used for systems which need Jacobian and residual evaluations, such as nonlinear systems being solved by Newton or quasi-Newton methods.
Similarly to MooseVariableFVReal, this variable describes a field which has been discretized using the cell-centered finite volume method and can be evaluated using the functor system in MOOSE.
Example Input File Syntax
To create a MooseLinearVariableFVReal
, users can do the following in their input files:
[Variables]
[v]
type = MooseLinearVariableFVReal
[]
[]
Note that unlike the case of finite element variables, the user needs to explicitly define the type of the variable.
Input Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Unit:(no unit assumed)
Controllable:No
Description:Adds user-defined labels for accessing object parameters via control logic.
- eigenFalseTrue to make this variable an eigen variable
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:True to make this variable an eigen variable
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Set the enabled status of the MooseObject.
- outputsVector of output names where you would like to restrict the output of variables(s) associated with this object
C++ Type:std::vector<OutputName>
Unit:(no unit assumed)
Controllable:No
Description:Vector of output names where you would like to restrict the output of variables(s) associated with this object
- scalingSpecifies a scaling factor to apply to this variable
C++ Type:std::vector<double>
Unit:(no unit assumed)
Controllable:No
Description:Specifies a scaling factor to apply to this variable
Advanced Parameters
Input Files
- (test/tests/linearfvkernels/diffusion-reaction-advection/advection-diffusion-reaction-2d.i)
- (test/tests/linearfvkernels/block-restriction/block-restricted-adr.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/mms/linear-segregated/2d-vortex/2d-vortex.i)
- (test/tests/linearfvkernels/block-restriction/block-restricted-diffusion.i)
- (test/tests/linearfvkernels/block-restriction/block-restricted-diffusion-react.i)
- (test/tests/outputs/debug/show_execution_linear_fv_elemental.i)
- (test/tests/linearfvkernels/diffusion/diffusion-2d-rz.i)
- (test/tests/linearfvkernels/diffusion/diffusion-1d.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/lid-driven/linear-segregated/lid-driven-segregated.i)
- (test/tests/variables/linearfv/diffusion-1d-aux.i)
- (test/tests/transfers/multiapp_copy_transfer/linear_sys_to_aux/nonlinear_main.i)
- (test/tests/postprocessors/side_diffusive_flux_integral/side_diffusive_flux_integral_linear_fv.i)
- (test/tests/transfers/multiapp_copy_transfer/linear_sys_to_aux/linear_sub.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/2d/2d-velocity-pressure.i)
- (test/tests/linearfvkernels/anisotropic-diffusion/anisotropic-diffusion-2d.i)
- (test/tests/variables/linearfv/diffusion-1d-pp.i)
- (test/tests/linearfvkernels/diffusion-reaction-advection/advection-diffusion-reaction-1d.i)
- (test/tests/linearfvkernels/advection/advection-1d.i)
- (test/tests/linearfvkernels/advection/advection-2d-rz.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/3d/3d-velocity-pressure.i)
- (test/tests/linearfvkernels/diffusion/diffusion-2d.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/segregated-comparison/segregated-linear.i)
- (test/tests/linearfvkernels/reaction/reaction-1d.i)
- (test/tests/linearfvkernels/advection/advection-2d.i)
- (test/tests/outputs/debug/show_execution_linear_fv_flux.i)
(test/tests/linearfvkernels/diffusion-reaction-advection/advection-diffusion-reaction-2d.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 2
nx = 2
ny = 1
ymax = 0.5
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[LinearFVKernels]
[diffusion]
type = LinearFVDiffusion
variable = u
diffusion_coeff = diff_coeff_func
use_nonorthogonal_correction = false
[]
[advection]
type = LinearFVAdvection
variable = u
velocity = "0.5 0 0"
advected_interp_method = average
[]
[reaction]
type = LinearFVReaction
variable = u
coeff = coeff_func
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
inactive = "outflow"
[dir]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "left right top bottom"
functor = analytic_solution
[]
[outflow]
type = LinearFVAdvectionDiffusionOutflowBC
variable = u
boundary = "right"
use_two_term_expansion = true
[]
[]
[Functions]
[diff_coeff_func]
type = ParsedFunction
expression = '1.0+0.5*x*y'
[]
[coeff_func]
type = ParsedFunction
expression = '1.0+1.0/(1+x*y)'
[]
[source_func]
type = ParsedFunction
expression = '-1.0*x*pi*sin((1/2)*x*pi)*cos(2*y*pi) - 0.25*y*pi*sin(2*y*pi)*cos((1/2)*x*pi) + (1.0 + 1.0/(x*y + 1))*(sin((1/2)*x*pi)*sin(2*y*pi) + 1.5) + (17/4)*pi^2*(0.5*x*y + 1.0)*sin((1/2)*x*pi)*sin(2*y*pi) + 0.25*pi*sin(2*y*pi)*cos((1/2)*x*pi)'
[]
[analytic_solution]
type = ParsedFunction
expression = 'sin((1/2)*x*pi)*sin(2*y*pi) + 1.5'
[]
[]
[Postprocessors]
[h]
type = AverageElementSize
execute_on = FINAL
[]
[error]
type = ElementL2FunctorError
approximate = u
exact = analytic_solution
execute_on = FINAL
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
[csv]
type = CSV
execute_on = FINAL
[]
[]
(test/tests/linearfvkernels/block-restriction/block-restricted-adr.i)
[Mesh]
[cmg]
type = CartesianMeshGenerator
dim = 2
dx = '0.1 1 0.1'
dy = '0.1 0.5 0.1'
ix = '1 2 1'
iy = '1 1 1'
subdomain_id = '1 1 1 1 2 3 1 1 1'
[]
[transform]
type = TransformGenerator
input = cmg
transform = TRANSLATE
vector_value = '-0.1 -0.1 0.0'
[]
[create_sides]
type = SideSetsBetweenSubdomainsGenerator
input = transform
new_boundary = sides
primary_block = 2
paired_block = 1
[]
[create_outlet]
type = SideSetsBetweenSubdomainsGenerator
input = create_sides
new_boundary = outlet
primary_block = 2
paired_block = 3
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
block = 2
[]
[]
[LinearFVKernels]
[diffusion]
type = LinearFVDiffusion
variable = u
diffusion_coeff = diff_coeff_func
use_nonorthogonal_correction = false
[]
[advection]
type = LinearFVAdvection
variable = u
velocity = "0.5 0 0"
advected_interp_method = average
[]
[reaction]
type = LinearFVReaction
variable = u
coeff = coeff_func
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
inactive = "outflow"
[dir]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "sides outlet"
functor = analytic_solution
[]
[outflow]
type = LinearFVAdvectionDiffusionOutflowBC
variable = u
boundary = "right"
use_two_term_expansion = true
[]
[]
[Functions]
[diff_coeff_func]
type = ParsedFunction
expression = '1.0+0.5*x*y'
[]
[coeff_func]
type = ParsedFunction
expression = '1.0+1.0/(1+x*y)'
[]
[source_func]
type = ParsedFunction
expression = '-1.0*x*pi*sin((1/2)*x*pi)*cos(2*y*pi) - 0.25*y*pi*sin(2*y*pi)*cos((1/2)*x*pi) + (1.0 + 1.0/(x*y + 1))*(sin((1/2)*x*pi)*sin(2*y*pi) + 1.5) + (17/4)*pi^2*(0.5*x*y + 1.0)*sin((1/2)*x*pi)*sin(2*y*pi) + 0.25*pi*sin(2*y*pi)*cos((1/2)*x*pi)'
[]
[analytic_solution]
type = ParsedFunction
expression = 'sin((1/2)*x*pi)*sin(2*y*pi) + 1.5'
[]
[]
[Postprocessors]
[h]
type = AverageElementSize
execute_on = FINAL
block = 2
[]
[error]
type = ElementL2FunctorError
approximate = u
exact = analytic_solution
execute_on = FINAL
block = 2
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
[csv]
type = CSV
execute_on = FINAL
[]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/mms/linear-segregated/2d-vortex/2d-vortex.i)
mu = 1
rho = 1
advected_interp_method = 'average'
[Problem]
linear_sys_names = 'u_system v_system pressure_system'
previous_nl_solution_required = true
[]
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 2
nx = 2
ny = 2
[]
[]
[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.0
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
[]
[]
[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_forcing]
type = LinearFVSource
variable = vel_x
source_density = forcing_u
[]
[v_forcing]
type = LinearFVSource
variable = vel_y
source_density = forcing_v
[]
[p_diffusion]
type = LinearFVAnisotropicDiffusion
variable = pressure
diffusion_tensor = Ainv
use_nonorthogonal_correction = false
use_nonorthogonal_correction_on_boundary = false
[]
[HbyA_divergence]
type = LinearFVDivergence
variable = pressure
face_flux = HbyA
force_boundary_execution = true
[]
[]
[LinearFVBCs]
[no-slip-wall-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left right top bottom'
variable = vel_x
functor = '0'
[]
[no-slip-wall-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left right top bottom'
variable = vel_y
functor = '0'
[]
[pressure-extrapolation]
type = LinearFVExtrapolatedPressureBC
boundary = 'left right top bottom'
variable = pressure
use_two_term_expansion = true
[]
[]
[Functions]
[exact_u]
type = ParsedFunction
expression = 'x^2*(1-x)^2*(2*y-6*y^2+4*y^3)'
[]
[exact_v]
type = ParsedFunction
expression = '-y^2*(1-y)^2*(2*x-6*x^2+4*x^3)'
[]
[exact_p]
type = ParsedFunction
expression = 'x*(1-x)'
[]
[forcing_u]
type = ParsedFunction
expression = '-4*mu*(-1+2*y)*(y^2-6*x*y^2+6*x^2*y^2-y+6*x*y-6*x^2*y+3*x^2-6*x^3+3*x^4)+1-2*x+rho*4*x^3'
'*y^2*(2*y^2-2*y+1)*(y-1)^2*(-1+2*x)*(x-1)^3'
symbol_names = 'mu rho'
symbol_values = '${mu} ${rho}'
[]
[forcing_v]
type = ParsedFunction
expression = '4*mu*(-1+2*x)*(x^2-6*y*x^2+6*x^2*y^2-x+6*x*y-6*x*y^2+3*y^2-6*y^3+3*y^4)+rho*4*y^3*x^2*(2'
'*x^2-2*x+1)*(x-1)^2*(-1+2*y)*(y-1)^3'
symbol_names = 'mu rho'
symbol_values = '${mu} ${rho}'
[]
[]
[Executioner]
type = SIMPLE
momentum_l_abs_tol = 1e-8
pressure_l_abs_tol = 1e-8
momentum_l_tol = 0
pressure_l_tol = 0
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
momentum_equation_relaxation = 0.8
pressure_variable_relaxation = 0.3
num_iterations = 2000
pressure_absolute_tolerance = 1e-8
momentum_absolute_tolerance = 1e-8
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'
print_fields = false
pin_pressure = true
pressure_pin_value = 0.25
pressure_pin_point = '0.5 0.5 0.0'
[]
[Outputs]
exodus = true
[csv]
type = CSV
execute_on = FINAL
[]
[]
[Postprocessors]
[h]
type = AverageElementSize
outputs = 'csv'
execute_on = FINAL
[]
[L2u]
type = ElementL2FunctorError
approximate = vel_x
exact = exact_u
outputs = 'csv'
execute_on = FINAL
[]
[L2v]
type = ElementL2FunctorError
approximate = vel_y
exact = exact_v
outputs = 'csv'
execute_on = FINAL
[]
[L2p]
approximate = pressure
exact = exact_p
type = ElementL2FunctorError
outputs = 'csv'
execute_on = FINAL
[]
[]
(test/tests/linearfvkernels/block-restriction/block-restricted-diffusion.i)
[Mesh]
[cmg]
type = CartesianMeshGenerator
dim = 2
dx = '0.1 1 0.1'
dy = '0.1 0.5 0.1'
ix = '1 2 1'
iy = '1 1 1'
subdomain_id = '1 1 1 1 2 3 1 1 1'
[]
[transform]
type = TransformGenerator
input = cmg
transform = TRANSLATE
vector_value = '-0.1 -0.1 0.0'
[]
[create_sides]
type = SideSetsBetweenSubdomainsGenerator
input = transform
new_boundary = sides
primary_block = 2
paired_block = 1
[]
[create_outlet]
type = SideSetsBetweenSubdomainsGenerator
input = create_sides
new_boundary = outlet
primary_block = 2
paired_block = 3
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
block = 2
[]
[]
[LinearFVKernels]
[diffusion]
type = LinearFVDiffusion
variable = u
diffusion_coeff = diff_coeff_func
use_nonorthogonal_correction = false
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
[dir]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "sides outlet"
functor = analytic_solution
[]
[]
[Functions]
[diff_coeff_func]
type = ParsedFunction
expression = '1.0+0.5*x*y'
[]
[source_func]
type = ParsedFunction
expression = '-1.0*x*pi*sin(x*pi)*cos(2*y*pi) - 0.5*y*pi*sin(2*y*pi)*cos(x*pi) + 5*pi^2*(0.5*x*y + 1.0)*sin(x*pi)*sin(2*y*pi)'
[]
[analytic_solution]
type = ParsedFunction
expression = 'sin(x*pi)*sin(2*y*pi) + 1.5'
[]
[]
[Postprocessors]
[h]
type = AverageElementSize
execute_on = FINAL
block = 2
[]
[error]
type = ElementL2FunctorError
approximate = u
exact = analytic_solution
execute_on = FINAL
block = 2
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
[csv]
type = CSV
execute_on = FINAL
[]
[]
(test/tests/linearfvkernels/block-restriction/block-restricted-diffusion-react.i)
source=1
diff_coeff=2
reac_coeff=3
[Mesh]
[cmg]
type = CartesianMeshGenerator
dim = 1
dx = '0.5 0.5'
ix = '20 20'
subdomain_id = '1 2'
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[LinearFVKernels]
[diffusion]
type = LinearFVDiffusion
variable = u
diffusion_coeff = ${diff_coeff}
use_nonorthogonal_correction = false
block = 1
[]
[reaction]
type = LinearFVReaction
variable = u
coeff = ${reac_coeff}
block = 2
[]
[source]
type = LinearFVSource
variable = u
source_density = ${source}
[]
[]
[LinearFVBCs]
[dir]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "left"
functor = 0
[]
[]
[Functions]
[analytic_solution]
type = ParsedFunction
expression = 'if(x<0.5, -x*x*S/2/D+(S/C+0.5*0.5/2/D*S)/0.5*x, S/C)'
symbol_names = 'S D C'
symbol_values = '${source} ${diff_coeff} ${reac_coeff}'
[]
[]
[Postprocessors]
[h]
type = AverageElementSize
execute_on = FINAL
block = 2
[]
[error]
type = ElementL2FunctorError
approximate = u
exact = analytic_solution
execute_on = FINAL
block = 2
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
[csv]
type = CSV
execute_on = FINAL
[]
[]
(test/tests/outputs/debug/show_execution_linear_fv_elemental.i)
[Mesh]
[gen_mesh]
type = GeneratedMeshGenerator
dim = 1
xmin = 0
xmax = 10
nx = 50
[]
[left]
type = ParsedSubdomainMeshGenerator
input = 'gen_mesh'
combinatorial_geometry = 'x < 0.5'
block_id = '1'
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
[]
[]
[LinearFVKernels]
[reaction_1]
type = LinearFVReaction
variable = u
coeff = 1.5
block = 0
[]
[reaction_2]
type = LinearFVReaction
variable = u
coeff = 2.5
block = 1
[]
[source_1]
type = LinearFVSource
variable = u
source_density = 3.5
block = 0
[]
[source_2]
type = LinearFVSource
variable = u
source_density = 4.5
block = 1
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
[]
[Debug]
show_execution_order = 'NONLINEAR'
[]
(test/tests/linearfvkernels/diffusion/diffusion-2d-rz.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 2
nx = 2
ny = 1
ymax = 0.5
[]
coord_type = RZ
rz_coord_axis = Y
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[LinearFVKernels]
[diffusion]
type = LinearFVDiffusion
variable = u
diffusion_coeff = coeff_func
use_nonorthogonal_correction = true
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
[dir]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "left right top bottom"
functor = analytic_solution
[]
[]
[Functions]
[coeff_func]
type = ParsedFunction
expression = '1+0.5*x*y'
[]
[source_func]
type = ParsedFunction
expression = '-(-1.0*x^2*y*(1.5 - x^2) + x*(1.5 - x^2)*(-1.0*x*y - 2))/x - (-1.0*x^2*y*(1.5 - y^2) - 4*x*(1.5 - y^2)*(0.5*x*y + 1))/x'
[]
[analytic_solution]
type = ParsedFunction
expression = '(1.5-x*x)*(1.5-y*y)'
[]
[]
[Postprocessors]
[h]
type = AverageElementSize
execute_on = FINAL
[]
[error]
type = ElementL2FunctorError
approximate = u
exact = analytic_solution
execute_on = FINAL
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
[csv]
type = CSV
execute_on = FINAL
[]
[exo]
type = Exodus
execute_on = FINAL
[]
[]
(test/tests/linearfvkernels/diffusion/diffusion-1d.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 1
nx = 2
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[LinearFVKernels]
[diffusion]
type = LinearFVDiffusion
variable = u
diffusion_coeff = coeff_func
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
[dir]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "left right"
functor = analytic_solution
[]
[]
[Functions]
[coeff_func]
type = ParsedFunction
expression = '0.5*x'
[]
[source_func]
type = ParsedFunction
expression = '2*x'
[]
[analytic_solution]
type = ParsedFunction
expression = '1-x*x'
[]
[]
[Postprocessors]
[h]
type = AverageElementSize
execute_on = FINAL
[]
[error]
type = ElementL2FunctorError
approximate = u
exact = analytic_solution
execute_on = FINAL
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
[csv]
type = CSV
execute_on = FINAL
[]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/lid-driven/linear-segregated/lid-driven-segregated.i)
mu = .01
rho = 1
advected_interp_method = 'average'
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = .1
ymin = 0
ymax = .1
nx = 3
ny = 3
[]
[]
[Problem]
linear_sys_names = 'u_system v_system pressure_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.0
solver_sys = u_system
[]
[vel_y]
type = MooseLinearVariableFVReal
initial_condition = 0.0
solver_sys = v_system
[]
[pressure]
type = MooseLinearVariableFVReal
solver_sys = pressure_system
initial_condition = 0.2
[]
[]
[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'
[]
[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
[]
[]
[LinearFVBCs]
[top_x]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = vel_x
boundary = 'top'
functor = 1
[]
[no_slip_x]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = vel_x
boundary = 'left right bottom'
functor = 0
[]
[no_slip_y]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = vel_y
boundary = 'left right top bottom'
functor = 0
[]
[pressure-extrapolation]
type = LinearFVExtrapolatedPressureBC
boundary = 'left right top bottom'
variable = pressure
use_two_term_expansion = true
[]
[]
[Executioner]
type = SIMPLE
momentum_l_abs_tol = 1e-12
pressure_l_abs_tol = 1e-12
momentum_l_tol = 0
pressure_l_tol = 0
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
momentum_equation_relaxation = 0.8
pressure_variable_relaxation = 0.3
num_iterations = 100
pressure_absolute_tolerance = 1e-12
momentum_absolute_tolerance = 1e-12
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'
print_fields = false
pin_pressure = true
pressure_pin_value = 0.0
pressure_pin_point = '0.05 0.05 0.0'
[]
[Outputs]
exodus = true
csv = false
perf_graph = false
print_nonlinear_residuals = false
print_linear_residuals = true
[]
(test/tests/variables/linearfv/diffusion-1d-aux.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 1
nx = 10
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[AuxVariables]
[v_volume]
type = MooseLinearVariableFVReal
initial_condition = 50
[]
[v_functor]
type = MooseLinearVariableFVReal
initial_condition = 25
[]
[v_parsed]
type = MooseLinearVariableFVReal
initial_condition = 12.5
[]
[]
[AuxKernels]
[volume]
type = VolumeAux
variable = v_volume
[]
[functor]
type = FunctorAux
variable = v_functor
functor = u
[]
[parsed]
type = ParsedAux
variable = v_parsed
coupled_variables = 'v_volume v_functor'
expression = '0.5*v_volume+0.5*v_functor'
[]
[]
[LinearFVKernels]
[diffusion]
type = LinearFVDiffusion
variable = u
diffusion_coeff = coeff_func
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
[dir]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "left right"
functor = analytic_solution
[]
[]
[Functions]
[coeff_func]
type = ParsedFunction
expression = '0.5*x'
[]
[source_func]
type = ParsedFunction
expression = '2*x'
[]
[analytic_solution]
type = ParsedFunction
expression = '1-x*x'
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
execute_on = TIMESTEP_END
[]
(test/tests/transfers/multiapp_copy_transfer/linear_sys_to_aux/nonlinear_main.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 1
nx = 10
[]
[]
[Variables]
[u_main]
type = MooseVariableFVReal
[]
[]
[AuxVariables]
[transferred]
type = MooseLinearVariableFVReal
[]
[]
[Transfers]
[copy]
type = MultiAppCopyTransfer
from_multi_app = linear_sub
source_variable = u
variable = transferred
[]
[]
[MultiApps]
[linear_sub]
type = FullSolveMultiApp
input_files = 'linear_sub.i'
[]
[]
[FVKernels]
[diff]
type = FVDiffusion
variable = u_main
coeff = 2
[]
[]
[FVBCs]
[left]
type = FVDirichletBC
variable = u_main
boundary = left
value = 0
[]
[./right]
type = FVDirichletBC
variable = u_main
boundary = right
value = 1
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/postprocessors/side_diffusive_flux_integral/side_diffusive_flux_integral_linear_fv.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
solver_sys = 'u_sys'
type = MooseLinearVariableFVReal
[]
[]
[LinearFVKernels]
[diff]
type = LinearFVDiffusion
variable = u
diffusion_coeff = diffusivity
[]
[]
[LinearFVBCs]
[left]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = left
functor = 0
[]
[right]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = right
functor = 1
[]
[]
[FunctorMaterials]
[mat_props]
type = GenericFunctorMaterial
prop_names = diffusivity
prop_values = 1
[]
[]
[AuxVariables]
[grad]
family = MONOMIAL_VEC
order = CONSTANT
[]
[]
[AuxKernels]
[grad]
variable = grad
functor = u
type = FunctorElementalGradientAux
[]
[]
[Postprocessors]
[avg_left]
type = SideAverageValue
variable = u
boundary = left
[]
[avg_right]
type = SideAverageValue
variable = u
boundary = right
[]
[avg_flux_left]
type = SideDiffusiveFluxAverage
variable = u
boundary = left
functor_diffusivity = diffusivity
[]
[avg_flux_right]
type = SideDiffusiveFluxAverage
variable = u
boundary = right
functor_diffusivity = diffusivity
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(test/tests/transfers/multiapp_copy_transfer/linear_sys_to_aux/linear_sub.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 1
nx = 10
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[LinearFVKernels]
[diffusion]
type = LinearFVDiffusion
variable = u
diffusion_coeff = coeff_func
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
[dir]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "left right"
functor = analytic_solution
[]
[]
[Functions]
[coeff_func]
type = ParsedFunction
expression = '0.5*x'
[]
[source_func]
type = ParsedFunction
expression = '2*x'
[]
[analytic_solution]
type = ParsedFunction
expression = '1-x*x'
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
execute_on = TIMESTEP_END
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/2d/2d-velocity-pressure.i)
mu = 2.6
rho = 1.0
advected_interp_method = 'average'
[Mesh]
[mesh]
type = CartesianMeshGenerator
dim = 2
dx = '0.3'
dy = '0.3'
ix = '3'
iy = '3'
[]
[]
[Problem]
linear_sys_names = 'u_system v_system pressure_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
[]
[]
[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'
[]
[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
[]
[]
[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
[]
[]
[Executioner]
type = SIMPLE
momentum_l_abs_tol = 1e-10
pressure_l_abs_tol = 1e-10
momentum_l_tol = 0
pressure_l_tol = 0
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
momentum_equation_relaxation = 0.8
pressure_variable_relaxation = 0.3
num_iterations = 50
pressure_absolute_tolerance = 1e-10
momentum_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'
print_fields = false
[]
[Outputs]
exodus = true
[]
(test/tests/linearfvkernels/anisotropic-diffusion/anisotropic-diffusion-2d.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 2
nx = 2
ny = 1
ymax = 0.5
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[LinearFVKernels]
[diffusion]
type = LinearFVAnisotropicDiffusion
variable = u
diffusion_tensor = diffusivity_tensor
use_nonorthogonal_correction = false
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
[dir]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "left right top bottom"
functor = analytic_solution
[]
[]
[FunctorMaterials]
[diff_tensor]
type = GenericVectorFunctorMaterial
prop_names = diffusivity_tensor
prop_values = 'coeff_func_x coeff_func_y 0.0'
[]
[]
[Functions]
[coeff_func_x]
type = ParsedFunction
expression = '1+0.5*x*y'
[]
[coeff_func_y]
type = ParsedFunction
expression = '1+x*y'
[]
[source_func]
type = ParsedFunction
expression = '(1.5-y*y)*(2+2*x*y)+(1.5-x*x)*(2+4*x*y)'
[]
[analytic_solution]
type = ParsedFunction
expression = '(1.5-x*x)*(1.5-y*y)'
[]
[]
[Postprocessors]
[h]
type = AverageElementSize
execute_on = FINAL
[]
[error]
type = ElementL2FunctorError
approximate = u
exact = analytic_solution
execute_on = FINAL
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
number_of_iterations = 1
[]
[Outputs]
[csv]
type = CSV
execute_on = FINAL
[]
[]
(test/tests/variables/linearfv/diffusion-1d-pp.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 1
nx = 50
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[LinearFVKernels]
[diffusion]
type = LinearFVDiffusion
variable = u
diffusion_coeff = coeff_func
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
[dir]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "left right"
functor = analytic_solution
[]
[]
[Functions]
[coeff_func]
type = ParsedFunction
expression = '0.5*x'
[]
[source_func]
type = ParsedFunction
expression = '2*x'
[]
[analytic_solution]
type = ParsedFunction
expression = '1-x*x'
[]
[]
[Postprocessors]
[average]
type = ElementAverageValue
variable = u
execute_on = FINAL
outputs = csv
[]
[min]
type = ElementExtremeValue
variable = u
value_type = min
execute_on = FINAL
outputs = csv
[]
[max]
type = ElementExtremeValue
variable = u
value_type = max
execute_on = FINAL
outputs = csv
[]
[num_dofs]
type = NumDOFs
execute_on = FINAL
outputs = csv
[]
[elem_value]
type = ElementalVariableValue
variable = u
elementid = 10
execute_on = FINAL
outputs = csv
[]
[point_value]
type = PointValue
variable = u
point = '0.33333 0 0'
execute_on = FINAL
outputs = csv
[]
[]
[VectorPostprocessors]
[line-sample]
type = LineValueSampler
variable = u
start_point = '0.13333 0 0'
end_point = '0.766666 0 0'
num_points = 9
sort_by = x
execute_on = FINAL
outputs = vpp_csv
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
[csv]
type = CSV
execute_on = FINAL
[]
[vpp_csv]
type = CSV
execute_on = FINAL
[]
[]
(test/tests/linearfvkernels/diffusion-reaction-advection/advection-diffusion-reaction-1d.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 1
nx = 2
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[LinearFVKernels]
[diffusion]
type = LinearFVDiffusion
variable = u
diffusion_coeff = diff_coeff_func
use_nonorthogonal_correction = false
[]
[advection]
type = LinearFVAdvection
variable = u
velocity = "0.5 0 0"
advected_interp_method = average
[]
[reaction]
type = LinearFVReaction
variable = u
coeff = coeff_func
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
inactive = "outflow"
[dir]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "left right"
functor = analytic_solution
[]
[outflow]
type = LinearFVAdvectionDiffusionOutflowBC
variable = u
boundary = "right"
use_two_term_expansion = true
[]
[]
[Functions]
[diff_coeff_func]
type = ParsedFunction
expression = '1+0.5*x'
[]
[coeff_func]
type = ParsedFunction
expression = '1+1/(1+x)'
[]
[source_func]
type = ParsedFunction
expression = '(1+1/(x+1))*(sin(pi/2*x)+1.5)+0.25*pi*pi*(0.5*x+1)*sin(pi/2*x)'
[]
[analytic_solution]
type = ParsedFunction
expression = 'sin(pi/2*x)+1.5'
[]
[]
[Postprocessors]
[h]
type = AverageElementSize
execute_on = FINAL
[]
[error]
type = ElementL2FunctorError
approximate = u
exact = analytic_solution
execute_on = FINAL
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
[csv]
type = CSV
execute_on = FINAL
[]
[]
(test/tests/linearfvkernels/advection/advection-1d.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 1
nx = 2
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[LinearFVKernels]
[advection]
type = LinearFVAdvection
variable = u
velocity = "0.5 0 0"
advected_interp_method = upwind
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
[inflow]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "left"
functor = analytic_solution
[]
[outflow]
type = LinearFVAdvectionDiffusionOutflowBC
variable = u
boundary = "right"
use_two_term_expansion = false
[]
[]
[Functions]
[source_func]
type = ParsedFunction
expression = '0.5*x'
[]
[analytic_solution]
type = ParsedFunction
expression = '0.5+0.5*x*x'
[]
[]
[Postprocessors]
[error]
type = ElementL2FunctorError
approximate = u
exact = analytic_solution
execute_on = FINAL
[]
[h]
type = AverageElementSize
execute_on = FINAL
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
number_of_iterations = 1
[]
[Outputs]
[csv]
type = CSV
execute_on = FINAL
[]
[]
(test/tests/linearfvkernels/advection/advection-2d-rz.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 2
nx = 2
ny= 1
ymax = 0.5
[]
coord_type = RZ
rz_coord_axis = Y
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[LinearFVKernels]
[advection]
type = LinearFVAdvection
variable = u
velocity = "0.0 0.5 0"
advected_interp_method = average
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
[inflow]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "left right bottom"
functor = analytic_solution
[]
[outflow]
type = LinearFVAdvectionDiffusionOutflowBC
variable = u
boundary = "top"
use_two_term_expansion = true
[]
[]
[Functions]
[source_func]
type = ParsedFunction
expression = '1.0*pi*sin(x*pi)*cos(2*y*pi)'
[]
[analytic_solution]
type = ParsedFunction
expression = 'sin(x*pi)*sin(2*y*pi) + 1.5'
[]
[]
[Postprocessors]
[error]
type = ElementL2FunctorError
approximate = u
exact = analytic_solution
execute_on = FINAL
[]
[h]
type = AverageElementSize
execute_on = FINAL
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
number_of_iterations = 2
petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_shift_amount'
petsc_options_value = 'lu NONZERO 1e-10'
[]
[Outputs]
[csv]
type = CSV
execute_on = FINAL
[]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/3d/3d-velocity-pressure.i)
mu = 2.6
rho = 1.0
advected_interp_method = 'average'
[Mesh]
[mesh]
type = CartesianMeshGenerator
dim = 3
dx = '0.3'
dy = '0.3'
dz = '0.3'
ix = '3'
iy = '3'
iz = '3'
[]
[]
[Problem]
linear_sys_names = 'u_system v_system w_system pressure_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = RhieChowMassFlux
u = vel_x
v = vel_y
w = vel_z
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
[]
[vel_z]
type = MooseLinearVariableFVReal
solver_sys = w_system
initial_condition = 0.0
[]
[pressure]
type = MooseLinearVariableFVReal
solver_sys = pressure_system
initial_condition = 0.2
[]
[]
[LinearFVKernels]
[u_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_x
advected_interp_method = ${advected_interp_method}
mu = ${mu}
u = vel_x
v = vel_y
w = vel_z
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
w = vel_z
momentum_component = 'y'
rhie_chow_user_object = 'rc'
use_nonorthogonal_correction = false
[]
[w_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_z
advected_interp_method = ${advected_interp_method}
mu = ${mu}
u = vel_x
v = vel_y
w = vel_z
momentum_component = 'z'
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'
[]
[w_pressure]
type = LinearFVMomentumPressure
variable = vel_z
pressure = pressure
momentum_component = 'z'
[]
[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
[]
[]
[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'
[]
[inlet-w]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_z
functor = '0.0'
[]
[walls-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom back front'
variable = vel_x
functor = 0.0
[]
[walls-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom back front'
variable = vel_y
functor = 0.0
[]
[walls-w]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom back front'
variable = vel_z
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
[]
[outlet_w]
type = LinearFVAdvectionDiffusionOutflowBC
variable = vel_z
use_two_term_expansion = false
boundary = right
[]
[]
[Executioner]
type = SIMPLE
momentum_l_abs_tol = 1e-10
pressure_l_abs_tol = 1e-10
momentum_l_tol = 0
pressure_l_tol = 0
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system w_system'
pressure_system = 'pressure_system'
momentum_equation_relaxation = 0.8
pressure_variable_relaxation = 0.3
num_iterations = 100
pressure_absolute_tolerance = 1e-10
momentum_absolute_tolerance = 1e-10
momentum_petsc_options_iname = '-pc_type -pc_hypre_type -pc_hypre_boomeramg_agg_nl -pc_hypre_boomeramg_agg_num_paths -pc_hypre_boomeramg_truncfactor -pc_hypre_boomeramg_strong_threshold -pc_hypre_boomeramg_coarsen_type -pc_hypre_boomeramg_interp_type'
momentum_petsc_options_value = 'hypre boomeramg 4 1 0.1 0.6 HMIS ext+i'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type -pc_hypre_boomeramg_agg_nl -pc_hypre_boomeramg_agg_num_paths -pc_hypre_boomeramg_truncfactor -pc_hypre_boomeramg_strong_threshold -pc_hypre_boomeramg_coarsen_type -pc_hypre_boomeramg_interp_type'
pressure_petsc_options_value = 'hypre boomeramg 2 1 0.1 0.6 HMIS ext+i'
print_fields = false
[]
[Outputs]
exodus = true
[]
(test/tests/linearfvkernels/diffusion/diffusion-2d.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 2
nx = 2
ny = 1
ymax = 0.5
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[LinearFVKernels]
[diffusion]
type = LinearFVDiffusion
variable = u
diffusion_coeff = coeff_func
use_nonorthogonal_correction = false
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
[dir]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "left right top bottom"
functor = analytic_solution
[]
[]
[Functions]
[coeff_func]
type = ParsedFunction
expression = '1+0.5*x*y'
[]
[source_func]
type = ParsedFunction
expression = '2*(1.5-y*y)+2*x*y*(1.5-y*y)+2*(1.5-x*x)+2*x*y*(1.5-x*x)'
[]
[analytic_solution]
type = ParsedFunction
expression = '(1.5-x*x)*(1.5-y*y)'
[]
[]
[Postprocessors]
[h]
type = AverageElementSize
execute_on = FINAL
[]
[error]
type = ElementL2FunctorError
approximate = u
exact = analytic_solution
execute_on = FINAL
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
[csv]
type = CSV
execute_on = FINAL
[]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/segregated-comparison/segregated-linear.i)
mu = 2.6
rho = 1.0
advected_interp_method = 'average'
[Mesh]
[mesh]
type = CartesianMeshGenerator
dim = 2
dx = '0.3'
dy = '0.3'
ix = '3'
iy = '3'
[]
[]
[Problem]
linear_sys_names = 'u_system v_system pressure_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
[]
[]
[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'
[]
[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
[]
[]
[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
[]
[]
[Executioner]
type = SIMPLE
momentum_l_abs_tol = 1e-10
pressure_l_abs_tol = 1e-10
momentum_l_tol = 0
pressure_l_tol = 0
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
momentum_equation_relaxation = 0.8
pressure_variable_relaxation = 0.3
num_iterations = 2
pressure_absolute_tolerance = 1e-10
momentum_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'
print_fields = false
[]
[Outputs]
exodus = true
execute_on = FINAL
[]
(test/tests/linearfvkernels/reaction/reaction-1d.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 1
nx = 10
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[LinearFVKernels]
[reaction]
type = LinearFVReaction
variable = u
coeff = coeff_func
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[Functions]
[coeff_func]
type = ParsedFunction
expression = '1+sin(x)'
[]
[source_func]
type = ParsedFunction
expression = '(1+sin(x))*(1+cos(x))'
[]
[analytic_solution]
type = ParsedFunction
expression = '1+cos(x)'
[]
[]
[Postprocessors]
[l2error]
type = ElementL2FunctorError
approximate = u
exact = analytic_solution
execute_on = FINAL
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
[]
[Outputs]
[exodus]
type = Exodus
execute_on = FINAL
[]
[]
(test/tests/linearfvkernels/advection/advection-2d.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 2
nx = 2
ny= 1
ymax = 0.5
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[LinearFVKernels]
[advection]
type = LinearFVAdvection
variable = u
velocity = "0.5 0 0"
advected_interp_method = upwind
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
[inflow]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "left top bottom"
functor = analytic_solution
[]
[outflow]
type = LinearFVAdvectionDiffusionOutflowBC
variable = u
boundary = "right"
use_two_term_expansion = false
[]
[]
[Functions]
[source_func]
type = ParsedFunction
expression = '0.5*pi*sin(2*y*pi)*cos(x*pi)'
[]
[analytic_solution]
type = ParsedFunction
expression = 'sin(x*pi)*sin(2*y*pi) + 1.5'
[]
[]
[Postprocessors]
[error]
type = ElementL2FunctorError
approximate = u
exact = analytic_solution
execute_on = FINAL
[]
[h]
type = AverageElementSize
execute_on = FINAL
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
number_of_iterations = 1
petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_shift_amount'
petsc_options_value = 'lu NONZERO 1e-10'
[]
[Outputs]
[csv]
type = CSV
execute_on = FINAL
[]
[]
(test/tests/outputs/debug/show_execution_linear_fv_flux.i)
[Mesh]
[gen_mesh]
type = GeneratedMeshGenerator
dim = 1
xmin = 0
xmax = 10
nx = 50
[]
[left]
type = ParsedSubdomainMeshGenerator
input = 'gen_mesh'
combinatorial_geometry = 'x < 0.5'
block_id = '1'
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
[]
[]
[LinearFVKernels]
[diffusion_1]
type = LinearFVDiffusion
variable = u
diffusion_coeff = 1.5
block = 0
[]
[diffusion_2]
type = LinearFVDiffusion
variable = u
diffusion_coeff = 2.5
block = 1
[]
[source_1]
type = LinearFVSource
variable = u
source_density = 3.5
block = 0
[]
[source_2]
type = LinearFVSource
variable = u
source_density = 4.5
block = 1
[]
[]
[LinearFVBCs]
[left_bc]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = 'left right'
functor = 0
[]
[]
[Executioner]
type = LinearPicardSteady
linear_systems_to_solve = u_sys
[]
[Debug]
show_execution_order = 'NONLINEAR'
[]