- boundaryThe list of boundary IDs from the mesh where this object applies
C++ Type:std::vector<BoundaryName>
Unit:(no unit assumed)
Controllable:No
Description:The list of boundary IDs from the mesh where this object applies
- eqnThe equation you're solving.
C++ Type:MooseEnum
Unit:(no unit assumed)
Controllable:No
Description:The equation you're solving.
- fpFluid properties userobject
C++ Type:UserObjectName
Unit:(no unit assumed)
Controllable:No
Description:Fluid properties userobject
- variableThe name of the variable that this boundary condition applies to
C++ Type:NonlinearVariableName
Unit:(no unit assumed)
Controllable:No
Description:The name of the variable that this boundary condition applies to
PCNSFVStrongBC
Computes the residual of advective term using finite volume method.
Overview
This object accepts functions describing boundary values for pressure, fluid temperature, and superficial velocity. The superficial velocity functions can also be used to supply superficial momentum information instead by setting velocity_function_includes_rho = true
. If no function is provided for a quantity, the boundary value of that quantity will be determined by extrapolating from the neighboring boundary cell centroid using cell centroid value and gradient information. From this mix of user-provided functions and extrapolated boundary values for pressure, fluid temperature, and superficial velocity/momentum, the fluxes for mass, momentum, energy, and even passive scalars can be computed.
Input Parameters
- T_fluidAn optional name of a function for the fluid temperature. If not provided then the fluid temperature will be treated implicitly (e.g. we will use the interior value
C++ Type:FunctionName
Unit:(no unit assumed)
Controllable:No
Description:An optional name of a function for the fluid temperature. If not provided then the fluid temperature will be treated implicitly (e.g. we will use the interior value
- displacementsThe displacements
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:The displacements
- mfr_postprocessorA postprocessor used for outputting mass flow rates on the same boundary this object acts on
C++ Type:UserObjectName
Unit:(no unit assumed)
Controllable:No
Description:A postprocessor used for outputting mass flow rates on the same boundary this object acts on
- momentum_componentThe component of the momentum equation that this kernel applies to.
C++ Type:MooseEnum
Unit:(no unit assumed)
Options:x, y, z
Controllable:No
Description:The component of the momentum equation that this kernel applies to.
- pressureAn optional name of a function for the pressure. If not provided then the pressure will be treated implicitly (e.g. we will use the interior value
C++ Type:FunctionName
Unit:(no unit assumed)
Controllable:No
Description:An optional name of a function for the pressure. If not provided then the pressure will be treated implicitly (e.g. we will use the interior value
- prop_getter_suffixAn optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:An optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
- scalarA function describing the value of the scalar at the boundary. If this function is not provided, then the interior value will be used.
C++ Type:FunctionName
Unit:(no unit assumed)
Controllable:No
Description:A function describing the value of the scalar at the boundary. If this function is not provided, then the interior value will be used.
- superficial_velocityAn optional name of a vector function for the velocity. If not provided then the superficial velocity will be treated implicitly (e.g. we will use the interior value
C++ Type:FunctionName
Unit:(no unit assumed)
Controllable:No
Description:An optional name of a vector function for the velocity. If not provided then the superficial velocity will be treated implicitly (e.g. we will use the interior value
- use_interpolated_stateFalseFor the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:For the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
- velocity_function_includes_rhoFalseWhether the provided superficial velocity function actually includes multiplication by rho (e.g. the function is representative of momentum.
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Whether the provided superficial velocity function actually includes multiplication by rho (e.g. the function is representative of momentum.
Optional Parameters
- absolute_value_vector_tagsThe tags for the vectors this residual object should fill with the absolute value of the residual contribution
C++ Type:std::vector<TagName>
Unit:(no unit assumed)
Controllable:No
Description:The tags for the vectors this residual object should fill with the absolute value of the residual contribution
- extra_matrix_tagsThe extra tags for the matrices this Kernel should fill
C++ Type:std::vector<TagName>
Unit:(no unit assumed)
Controllable:No
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<TagName>
Unit:(no unit assumed)
Controllable:No
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
Unit:(no unit assumed)
Options:nontime, system
Controllable:No
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
Unit:(no unit assumed)
Options:nontime, time
Controllable:No
Description:The tag for the vectors this Kernel should fill
Tagging 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.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Unit:(no unit assumed)
Controllable:Yes
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
Unit:(no unit assumed)
Controllable:No
Description:Determines whether this object is calculated using an implicit or explicit form
- 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
Unit:(no unit assumed)
Controllable:No
Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Advanced Parameters
Input Files
- (modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/varying-eps-basic-kt-primitive.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/heated-channel/transient-porous-kt-primitive.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/scalar_advection/mass-frac-advection.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/straight_channel_porosity_step/rotated-2d-bkt-function-porosity-mixed.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/straight_channel_porosity_step/dc.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/basic-conserved-pcnsfv-kt.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/straight_channel_porosity_step/rotated-2d-bkt-function-porosity.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/straight_channel_porosity_step/hllc.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/varying-eps-hllc.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/straight_channel_porosity_step/implicit-euler-basic-kt-primitive.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/varying-eps-basic-kt-mixed.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/basic-primitive-pcnsfv-kt.i)
(modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/varying-eps-basic-kt-primitive.i)
[GlobalParams]
fp = fp
limiter = 'central_difference'
two_term_boundary_expansion = true
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 1
xmin = .1
xmax = .6
nx = 2
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[Variables]
[pressure]
type = MooseVariableFVReal
[]
[sup_vel_x]
type = MooseVariableFVReal
[]
[T_fluid]
type = MooseVariableFVReal
[]
[]
[ICs]
[pressure]
type = FunctionIC
variable = pressure
function = 'exact_p'
[]
[sup_vel_x]
type = FunctionIC
variable = sup_vel_x
function = 'exact_sup_vel_x'
[]
[T_fluid]
type = FunctionIC
variable = T_fluid
function = 'exact_T'
[]
[]
[FVKernels]
[mass_advection]
type = PCNSFVKT
variable = pressure
eqn = "mass"
[]
[mass_fn]
type = FVBodyForce
variable = pressure
function = 'forcing_rho'
[]
[momentum_x_advection]
type = PCNSFVKT
variable = sup_vel_x
momentum_component = x
eqn = "momentum"
[]
[eps_grad]
type = PNSFVPGradEpsilon
variable = sup_vel_x
momentum_component = 'x'
epsilon_function = 'eps'
[]
[momentum_fn]
type = FVBodyForce
variable = sup_vel_x
function = 'forcing_rho_ud'
[]
[fluid_energy_advection]
type = PCNSFVKT
variable = T_fluid
eqn = "energy"
[]
[energy_fn]
type = FVBodyForce
variable = T_fluid
function = 'forcing_rho_et'
[]
[]
[FVBCs]
[mass_left]
variable = pressure
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'mass'
[]
[momentum_left]
variable = sup_vel_x
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'momentum'
momentum_component = 'x'
[]
[energy_left]
variable = T_fluid
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'energy'
[]
[mass_right]
variable = pressure
type = PCNSFVStrongBC
boundary = right
eqn = 'mass'
pressure = 'exact_p'
[]
[momentum_right]
variable = sup_vel_x
type = PCNSFVStrongBC
boundary = right
eqn = 'momentum'
momentum_component = 'x'
pressure = 'exact_p'
[]
[energy_right]
variable = T_fluid
type = PCNSFVStrongBC
boundary = right
eqn = 'energy'
pressure = 'exact_p'
[]
# help gradient reconstruction
[pressure_right]
type = FVFunctionDirichletBC
variable = pressure
function = exact_p
boundary = 'right'
[]
[sup_vel_x_left]
type = FVFunctionDirichletBC
variable = sup_vel_x
function = exact_sup_vel_x
boundary = 'left'
[]
[T_fluid_left]
type = FVFunctionDirichletBC
variable = T_fluid
function = exact_T
boundary = 'left'
[]
[]
[Materials]
[var_mat]
type = PorousPrimitiveVarMaterial
pressure = pressure
superficial_vel_x = sup_vel_x
T_fluid = T_fluid
porosity = porosity
[]
[porosity]
type = GenericFunctionMaterial
prop_names = 'porosity'
prop_values = 'eps'
[]
[]
[Functions]
[exact_rho]
type = ParsedFunction
expression = '3.48788261470924*cos(x)'
[]
[forcing_rho]
type = ParsedFunction
expression = '-3.83667087618017*sin(1.1*x)*cos(1.3*x) - 4.53424739912202*sin(1.3*x)*cos(1.1*x)'
[]
[exact_rho_ud]
type = ParsedFunction
expression = '3.48788261470924*cos(1.1*x)*cos(1.3*x)'
[]
[forcing_rho_ud]
type = ParsedFunction
expression = '(-(10.6975765229419*cos(1.5*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + (10.6975765229419*sin(x)*cos(1.5*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 16.0463647844128*sin(1.5*x)/cos(x))*cos(x))*cos(1.3*x) + 3.48788261470924*sin(x)*cos(1.1*x)^2*cos(1.3*x)/cos(x)^2 - 7.67334175236034*sin(1.1*x)*cos(1.1*x)*cos(1.3*x)/cos(x) - 4.53424739912202*sin(1.3*x)*cos(1.1*x)^2/cos(x)'
[]
[exact_rho_et]
type = ParsedFunction
expression = '26.7439413073546*cos(1.5*x)'
[]
[forcing_rho_et]
type = ParsedFunction
expression = '1.0*(3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.5*x))*sin(x)*cos(1.1*x)*cos(1.3*x)/cos(x)^2 - 1.1*(3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.5*x))*sin(1.1*x)*cos(1.3*x)/cos(x) - 1.3*(3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.5*x))*sin(1.3*x)*cos(1.1*x)/cos(x) + 1.0*(-(10.6975765229419*cos(1.5*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + (10.6975765229419*sin(x)*cos(1.5*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 16.0463647844128*sin(1.5*x)/cos(x))*cos(x) - 40.1159119610319*sin(1.5*x))*cos(1.1*x)*cos(1.3*x)/cos(x)'
[]
[exact_T]
type = ParsedFunction
expression = '0.0106975765229418*cos(1.5*x)/cos(x) - 0.000697576522941848*cos(1.1*x)^2/cos(x)^2'
[]
[exact_eps_p]
type = ParsedFunction
expression = '3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)*cos(1.3*x)'
[]
[exact_p]
type = ParsedFunction
expression = '3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)'
[]
[exact_sup_vel_x]
type = ParsedFunction
expression = '1.0*cos(1.1*x)*cos(1.3*x)/cos(x)'
[]
[eps]
type = ParsedFunction
expression = 'cos(1.3*x)'
[]
[exact_superficial_velocity]
type = ParsedVectorFunction
expression_x = '1.0*cos(1.1*x)*cos(1.3*x)/cos(x)'
[]
[]
[Executioner]
solve_type = NEWTON
type = Transient
num_steps = 1
dtmin = 1
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
nl_max_its = 50
line_search = bt
nl_rel_tol = 1e-12
nl_abs_tol = 1e-12
[]
[Outputs]
exodus = true
csv = true
[]
[Debug]
show_var_residual_norms = true
[]
[Postprocessors]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2pressure]
type = ElementL2Error
variable = pressure
function = exact_p
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2sup_vel_x]
variable = sup_vel_x
function = exact_sup_vel_x
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2T_fluid]
variable = T_fluid
function = exact_T
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/cns/heated-channel/transient-porous-kt-primitive.i)
p_initial=1.01e5
T=273.15
u_in=10
eps=1
superficial_vel_in=${fparse u_in * eps}
[GlobalParams]
fp = fp
limiter = 'vanLeer'
two_term_boundary_expansion = true
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 1
xmin = 0
xmax = 10
nx = 100
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[Variables]
[pressure]
type = MooseVariableFVReal
initial_condition = ${p_initial}
[]
[superficial_vel_x]
type = MooseVariableFVReal
initial_condition = ${superficial_vel_in}
[]
[temperature]
type = MooseVariableFVReal
initial_condition = ${T}
[]
[]
[AuxVariables]
[rho]
type = MooseVariableFVReal
[]
[superficial_rhou]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[rho]
type = ADMaterialRealAux
variable = rho
property = rho
execute_on = 'timestep_end'
[]
[superficial_rhou]
type = ADMaterialRealAux
variable = superficial_rhou
property = superficial_rhou
execute_on = 'timestep_end'
[]
[]
[FVKernels]
[mass_time]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rho_dt'
variable = pressure
[]
[mass_advection]
type = PCNSFVKT
variable = pressure
eqn = "mass"
[]
[momentum_time]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rhou_dt'
variable = superficial_vel_x
[]
[momentum_advection]
type = PCNSFVKT
variable = superficial_vel_x
eqn = "momentum"
momentum_component = 'x'
[]
[energy_time]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rho_et_dt'
variable = temperature
[]
[energy_advection]
type = PCNSFVKT
variable = temperature
eqn = "energy"
[]
[heat]
type = FVBodyForce
variable = temperature
value = 1e6
[]
[]
[FVBCs]
[rho_left]
type = PCNSFVStrongBC
boundary = 'left'
variable = pressure
superficial_velocity = 'superficial_vel_in'
T_fluid = ${T}
eqn = 'mass'
[]
[rhou_left]
type = PCNSFVStrongBC
boundary = 'left'
variable = superficial_vel_x
superficial_velocity = 'superficial_vel_in'
T_fluid = ${T}
eqn = 'momentum'
momentum_component = 'x'
[]
[rho_et_left]
type = PCNSFVStrongBC
boundary = 'left'
variable = temperature
superficial_velocity = 'superficial_vel_in'
T_fluid = ${T}
eqn = 'energy'
[]
[rho_right]
type = PCNSFVStrongBC
boundary = 'right'
variable = pressure
pressure = ${p_initial}
eqn = 'mass'
[]
[rhou_right]
type = PCNSFVStrongBC
boundary = 'right'
variable = superficial_vel_x
pressure = ${p_initial}
eqn = 'momentum'
momentum_component = 'x'
[]
[rho_et_right]
type = PCNSFVStrongBC
boundary = 'right'
variable = temperature
pressure = ${p_initial}
eqn = 'energy'
[]
# Use these to help create more accurate cell centered gradients for cells adjacent to boundaries
[T_left]
type = FVDirichletBC
variable = temperature
value = ${T}
boundary = 'left'
[]
[sup_vel_left]
type = FVDirichletBC
variable = superficial_vel_x
value = ${superficial_vel_in}
boundary = 'left'
[]
[p_right]
type = FVDirichletBC
variable = pressure
value = ${p_initial}
boundary = 'right'
[]
[]
[Functions]
[superficial_vel_in]
type = ParsedVectorFunction
expression_x = '${superficial_vel_in}'
[]
[]
[Materials]
[var_mat]
type = PorousPrimitiveVarMaterial
pressure = pressure
T_fluid = temperature
superficial_vel_x = superficial_vel_x
fp = fp
porosity = porosity
[]
[fluid_only]
type = GenericConstantMaterial
prop_names = 'porosity'
prop_values = '${eps}'
[]
[]
[Executioner]
solve_type = NEWTON
type = Transient
nl_max_its = 20
[TimeStepper]
type = IterationAdaptiveDT
dt = 5e-5
optimal_iterations = 10
[]
steady_state_detection = false
steady_state_tolerance = 1e-12
abort_on_solve_fail = false
end_time = 100
nl_abs_tol = 1e-8
dtmin = 5e-5
automatic_scaling = true
compute_scaling_once = false
verbose = true
petsc_options_iname = '-pc_type -pc_factor_mat_solver_type -pc_factor_shift_type -snes_linesearch_minlambda'
petsc_options_value = 'lu mumps NONZERO 1e-3 '
[]
[Outputs]
[exo]
type = Exodus
execute_on = 'final'
[]
[dof]
type = DOFMap
execute_on = 'initial'
[]
checkpoint = true
[]
[Debug]
show_var_residual_norms = true
[]
(modules/navier_stokes/test/tests/finite_volume/cns/scalar_advection/mass-frac-advection.i)
rho_initial=1.29
p_initial=1.01e5
T=273.15
gamma=1.4
e_initial=${fparse p_initial / (gamma - 1) / rho_initial}
et_initial=${e_initial}
rho_et_initial=${fparse rho_initial * et_initial}
v_in=1
[GlobalParams]
fp = fp
# retain behavior at time of test creation
two_term_boundary_expansion = false
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1
nx = 2
ymin = 0
ymax = 10
ny = 20
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[Variables]
[rho]
type = MooseVariableFVReal
initial_condition = ${rho_initial}
[]
[rho_u]
type = MooseVariableFVReal
initial_condition = 1e-15
[]
[rho_v]
type = MooseVariableFVReal
initial_condition = 1e-15
[]
[rho_et]
type = MooseVariableFVReal
initial_condition = ${rho_et_initial}
scaling = 1e-5
[]
[mass_frac]
type = MooseVariableFVReal
initial_condition = 1e-15
[]
[]
[AuxVariables]
[U_x]
type = MooseVariableFVReal
[]
[U_y]
type = MooseVariableFVReal
[]
[pressure]
type = MooseVariableFVReal
[]
[temperature]
type = MooseVariableFVReal
[]
[courant]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[U_x]
type = ADMaterialRealAux
variable = U_x
property = vel_x
execute_on = 'timestep_end'
[]
[U_y]
type = ADMaterialRealAux
variable = U_y
property = vel_y
execute_on = 'timestep_end'
[]
[pressure]
type = ADMaterialRealAux
variable = pressure
property = pressure
execute_on = 'timestep_end'
[]
[temperature]
type = ADMaterialRealAux
variable = temperature
property = T_fluid
execute_on = 'timestep_end'
[]
[courant]
type = Courant
variable = courant
u = U_x
v = U_y
[]
[]
[FVKernels]
[mass_time]
type = FVPorosityTimeDerivative
variable = rho
[]
[mass_advection]
type = PCNSFVKT
variable = rho
eqn = "mass"
[]
[momentum_time_x]
type = FVTimeKernel
variable = rho_u
[]
[momentum_advection_and_pressure_x]
type = PCNSFVKT
variable = rho_u
eqn = "momentum"
momentum_component = 'x'
[]
[momentum_time_y]
type = FVTimeKernel
variable = rho_v
[]
[momentum_advection_and_pressure_y]
type = PCNSFVKT
variable = rho_v
eqn = "momentum"
momentum_component = 'y'
[]
[energy_time]
type = FVPorosityTimeDerivative
variable = rho_et
[]
[energy_advection]
type = PCNSFVKT
variable = rho_et
eqn = "energy"
[]
[mass_frac_time]
type = PCNSFVDensityTimeDerivative
variable = mass_frac
rho = rho
[]
[mass_frac_advection]
type = PCNSFVKT
variable = mass_frac
eqn = "scalar"
[]
[]
[Functions]
[ud_in]
type = ParsedVectorFunction
expression_x = '0'
expression_y = '${v_in}'
[]
[]
[FVBCs]
[rho_bottom]
type = PCNSFVStrongBC
boundary = 'bottom'
variable = rho
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'mass'
[]
[rho_u_bottom]
type = PCNSFVStrongBC
boundary = 'bottom'
variable = rho_u
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'momentum'
momentum_component = 'x'
[]
[rho_v_bottom]
type = PCNSFVStrongBC
boundary = 'bottom'
variable = rho_v
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'momentum'
momentum_component = 'y'
[]
[rho_et_bottom]
type = PCNSFVStrongBC
boundary = 'bottom'
variable = rho_et
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'energy'
[]
[mass_frac_bottom]
type = PCNSFVStrongBC
boundary = 'bottom'
variable = mass_frac
superficial_velocity = 'ud_in'
T_fluid = ${T}
scalar = 1
eqn = 'scalar'
[]
[rho_top]
type = PCNSFVStrongBC
boundary = 'top'
variable = rho
pressure = ${p_initial}
eqn = 'mass'
[]
[rho_u_top]
type = PCNSFVStrongBC
boundary = 'top'
variable = rho_u
pressure = ${p_initial}
eqn = 'momentum'
momentum_component = 'x'
[]
[rho_v_top]
type = PCNSFVStrongBC
boundary = 'top'
variable = rho_v
pressure = ${p_initial}
eqn = 'momentum'
momentum_component = 'y'
[]
[rho_et_top]
type = PCNSFVStrongBC
boundary = 'top'
variable = rho_et
pressure = ${p_initial}
eqn = 'energy'
[]
[mass_frac_top]
type = PCNSFVStrongBC
boundary = 'top'
variable = mass_frac
pressure = ${p_initial}
eqn = 'scalar'
[]
[momentum_x_walls]
type = PCNSFVImplicitMomentumPressureBC
variable = rho_u
boundary = 'left right'
momentum_component = 'x'
[]
[momentum_y_walls]
type = PCNSFVImplicitMomentumPressureBC
variable = rho_v
boundary = 'left right'
momentum_component = 'y'
[]
[]
[Materials]
[var_mat]
type = PorousConservedVarMaterial
rho = rho
rho_et = rho_et
superficial_rhou = rho_u
superficial_rhov = rho_v
fp = fp
porosity = porosity
[]
[porosity]
type = GenericConstantMaterial
prop_names = 'porosity'
prop_values = '1'
[]
[]
[Executioner]
type = Transient
[TimeIntegrator]
type = ActuallyExplicitEuler
[]
steady_state_detection = true
steady_state_tolerance = 1e-12
abort_on_solve_fail = true
dt = 5e-4
num_steps = 25
[]
[Outputs]
[out]
type = Exodus
execute_on = 'initial timestep_end'
[]
[dof]
type = DOFMap
execute_on = 'initial'
[]
[]
[Debug]
show_var_residual_norms = true
[]
(modules/navier_stokes/test/tests/finite_volume/cns/straight_channel_porosity_step/rotated-2d-bkt-function-porosity-mixed.i)
p_initial=1.01e5
T=273.15
# u refers to the superficial velocity
u_in=1
rho_in=1.30524
sup_mom_y_in=${fparse u_in * rho_in}
user_limiter='upwind'
friction_coeff=10
[GlobalParams]
fp = fp
two_term_boundary_expansion = true
limiter = ${user_limiter}
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1
nx = 3
ymin = 0
ymax = 18
ny = 90
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[Variables]
[pressure]
type = MooseVariableFVReal
initial_condition = ${p_initial}
[]
[sup_mom_x]
type = MooseVariableFVReal
initial_condition = 1e-15
scaling = 1e-2
[]
[sup_mom_y]
type = MooseVariableFVReal
initial_condition = 1e-15
scaling = 1e-2
[]
[T_fluid]
type = MooseVariableFVReal
initial_condition = ${T}
scaling = 1e-5
[]
[]
[AuxVariables]
[vel_y]
type = MooseVariableFVReal
[]
[rho]
type = MooseVariableFVReal
[]
[eps]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[vel_y]
type = ADMaterialRealAux
variable = vel_y
property = vel_y
execute_on = 'timestep_end'
[]
[rho]
type = ADMaterialRealAux
variable = rho
property = rho
execute_on = 'timestep_end'
[]
[eps]
type = MaterialRealAux
variable = eps
property = porosity
execute_on = 'timestep_end'
[]
[]
[FVKernels]
[mass_time]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rho_dt'
variable = pressure
[]
[mass_advection]
type = PCNSFVKT
variable = pressure
eqn = "mass"
[]
[momentum_time]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rhou_dt'
variable = sup_mom_x
[]
[momentum_advection]
type = PCNSFVKT
variable = sup_mom_x
eqn = "momentum"
momentum_component = 'x'
[]
[eps_grad]
type = PNSFVPGradEpsilon
variable = sup_mom_x
momentum_component = 'x'
epsilon_function = 'eps'
[]
[drag]
type = PCNSFVMomentumFriction
variable = sup_mom_x
momentum_component = 'x'
Darcy_name = 'cl'
momentum_name = superficial_rhou
[]
[momentum_time_y]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rhov_dt'
variable = sup_mom_y
[]
[momentum_advection_y]
type = PCNSFVKT
variable = sup_mom_y
eqn = "momentum"
momentum_component = 'y'
[]
[eps_grad_y]
type = PNSFVPGradEpsilon
variable = sup_mom_y
momentum_component = 'y'
epsilon_function = 'eps'
[]
[drag_y]
type = PCNSFVMomentumFriction
variable = sup_mom_y
momentum_component = 'y'
Darcy_name = 'cl'
momentum_name = superficial_rhov
[]
[energy_time]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rho_et_dt'
variable = T_fluid
[]
[energy_advection]
type = PCNSFVKT
variable = T_fluid
eqn = "energy"
[]
[]
[FVBCs]
[rho_bottom]
type = PCNSFVStrongBC
boundary = 'bottom'
variable = pressure
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'mass'
velocity_function_includes_rho = true
[]
[rhou_bottom]
type = PCNSFVStrongBC
boundary = 'bottom'
variable = sup_mom_x
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'momentum'
momentum_component = 'x'
velocity_function_includes_rho = true
[]
[rhov_bottom]
type = PCNSFVStrongBC
boundary = 'bottom'
variable = sup_mom_y
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'momentum'
momentum_component = 'y'
velocity_function_includes_rho = true
[]
[rho_et_bottom]
type = PCNSFVStrongBC
boundary = 'bottom'
variable = T_fluid
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'energy'
velocity_function_includes_rho = true
[]
[rho_top]
type = PCNSFVStrongBC
boundary = 'top'
variable = pressure
pressure = ${p_initial}
eqn = 'mass'
[]
[rhou_top]
type = PCNSFVStrongBC
boundary = 'top'
variable = sup_mom_x
pressure = ${p_initial}
eqn = 'momentum'
momentum_component = 'x'
[]
[rhov_top]
type = PCNSFVStrongBC
boundary = 'top'
variable = sup_mom_y
pressure = ${p_initial}
eqn = 'momentum'
momentum_component = 'y'
[]
[rho_et_top]
type = PCNSFVStrongBC
boundary = 'top'
variable = T_fluid
pressure = ${p_initial}
eqn = 'energy'
[]
[wall_pressure_x]
type = PCNSFVImplicitMomentumPressureBC
momentum_component = 'x'
boundary = 'left right'
variable = sup_mom_x
[]
[wall_pressure_y]
type = PCNSFVImplicitMomentumPressureBC
momentum_component = 'y'
boundary = 'left right'
variable = sup_mom_y
[]
# Use these to help create more accurate cell centered gradients for cells adjacent to boundaries
[T_bottom]
type = FVDirichletBC
variable = T_fluid
value = ${T}
boundary = 'bottom'
[]
[sup_mom_x_bottom_and_walls]
type = FVDirichletBC
variable = sup_mom_x
value = 0
boundary = 'bottom left right'
[]
[sup_mom_y_walls]
type = FVDirichletBC
variable = sup_mom_y
value = 0
boundary = 'left right'
[]
[sup_mom_y_bottom]
type = FVDirichletBC
variable = sup_mom_y
value = ${sup_mom_y_in}
boundary = 'bottom'
[]
[p_top]
type = FVDirichletBC
variable = pressure
value = ${p_initial}
boundary = 'top'
[]
[]
[Functions]
[ud_in]
type = ParsedVectorFunction
expression_x = '0'
expression_y = '${sup_mom_y_in}'
[]
[eps]
type = ParsedFunction
expression = 'if(y < 2.8, 1,
if(y < 3.2, 1 - .5 / .4 * (y - 2.8),
if(y < 6.8, .5,
if(y < 7.2, .5 - .25 / .4 * (y - 6.8),
if(y < 10.8, .25,
if(y < 11.2, .25 + .25 / .4 * (y - 10.8),
if(y < 14.8, .5,
if(y < 15.2, .5 + .5 / .4 * (y - 14.8),
1))))))))'
[]
[]
[Materials]
[var_mat]
type = PorousMixedVarMaterial
pressure = pressure
T_fluid = T_fluid
superficial_rhou = sup_mom_x
superficial_rhov = sup_mom_y
fp = fp
porosity = porosity
[]
[porosity]
type = GenericFunctionMaterial
prop_names = 'porosity'
prop_values = 'eps'
[]
[ad_generic]
type = ADGenericConstantVectorMaterial
prop_names = 'cl'
prop_values = '${friction_coeff} ${friction_coeff} ${friction_coeff}'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
solve_type = NEWTON
line_search = 'bt'
type = Transient
nl_max_its = 20
[TimeStepper]
type = IterationAdaptiveDT
dt = 5e-5
optimal_iterations = 6
growth_factor = 1.2
[]
num_steps = 10000
end_time = 500
nl_abs_tol = 1e-7
petsc_options_iname = '-pc_type -pc_factor_mat_solver_type'
petsc_options_value = 'lu mumps'
[]
[Outputs]
[out]
type = Exodus
execute_on = 'final'
[]
checkpoint = true
[]
[Debug]
show_var_residual_norms = true
[]
(modules/navier_stokes/test/tests/finite_volume/cns/straight_channel_porosity_step/dc.i)
p_initial=1.01e5
T=273.15
# u refers to the superficial velocity
u_in=1
rho_in=1.30524
sup_mom_y_in=${fparse u_in * rho_in}
user_limiter='min_mod'
friction_coeff=10
[GlobalParams]
fp = fp
two_term_boundary_expansion = true
limiter = ${user_limiter}
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1
nx = 3
ymin = 0
ymax = 18
ny = 90
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[Variables]
[pressure]
type = MooseVariableFVReal
initial_condition = ${p_initial}
[]
[sup_mom_x]
type = MooseVariableFVReal
initial_condition = 1e-15
[]
[sup_mom_y]
type = MooseVariableFVReal
initial_condition = 1e-15
[]
[T_fluid]
type = MooseVariableFVReal
initial_condition = ${T}
[]
[]
[AuxVariables]
[vel_y]
type = MooseVariableFVReal
[]
[rho]
type = MooseVariableFVReal
[]
[eps]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[vel_y]
type = ADMaterialRealAux
variable = vel_y
property = vel_y
execute_on = 'timestep_end'
[]
[rho]
type = ADMaterialRealAux
variable = rho
property = rho
execute_on = 'timestep_end'
[]
[eps]
type = MaterialRealAux
variable = eps
property = porosity
execute_on = 'timestep_end'
[]
[]
[FVKernels]
[mass_time]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rho_dt'
variable = pressure
[]
[mass_advection]
type = PCNSFVKTDC
variable = pressure
eqn = "mass"
[]
[momentum_time]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rhou_dt'
variable = sup_mom_x
[]
[momentum_advection]
type = PCNSFVKTDC
variable = sup_mom_x
eqn = "momentum"
momentum_component = 'x'
[]
[eps_grad]
type = PNSFVPGradEpsilon
variable = sup_mom_x
momentum_component = 'x'
epsilon_function = 'eps'
[]
[drag]
type = PCNSFVMomentumFriction
variable = sup_mom_x
momentum_component = 'x'
Darcy_name = 'cl'
momentum_name = superficial_rhou
[]
[momentum_time_y]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rhov_dt'
variable = sup_mom_y
[]
[momentum_advection_y]
type = PCNSFVKTDC
variable = sup_mom_y
eqn = "momentum"
momentum_component = 'y'
[]
[eps_grad_y]
type = PNSFVPGradEpsilon
variable = sup_mom_y
momentum_component = 'y'
epsilon_function = 'eps'
[]
[drag_y]
type = PCNSFVMomentumFriction
variable = sup_mom_y
momentum_component = 'y'
Darcy_name = 'cl'
momentum_name = superficial_rhov
[]
[energy_time]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rho_et_dt'
variable = T_fluid
[]
[energy_advection]
type = PCNSFVKTDC
variable = T_fluid
eqn = "energy"
[]
[]
[FVBCs]
[rho_bottom]
type = PCNSFVStrongBC
boundary = 'bottom'
variable = pressure
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'mass'
velocity_function_includes_rho = true
[]
[rhou_bottom]
type = PCNSFVStrongBC
boundary = 'bottom'
variable = sup_mom_x
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'momentum'
momentum_component = 'x'
velocity_function_includes_rho = true
[]
[rhov_bottom]
type = PCNSFVStrongBC
boundary = 'bottom'
variable = sup_mom_y
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'momentum'
momentum_component = 'y'
velocity_function_includes_rho = true
[]
[rho_et_bottom]
type = PCNSFVStrongBC
boundary = 'bottom'
variable = T_fluid
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'energy'
velocity_function_includes_rho = true
[]
[rho_top]
type = PCNSFVStrongBC
boundary = 'top'
variable = pressure
pressure = ${p_initial}
eqn = 'mass'
[]
[rhou_top]
type = PCNSFVStrongBC
boundary = 'top'
variable = sup_mom_x
pressure = ${p_initial}
eqn = 'momentum'
momentum_component = 'x'
[]
[rhov_top]
type = PCNSFVStrongBC
boundary = 'top'
variable = sup_mom_y
pressure = ${p_initial}
eqn = 'momentum'
momentum_component = 'y'
[]
[rho_et_top]
type = PCNSFVStrongBC
boundary = 'top'
variable = T_fluid
pressure = ${p_initial}
eqn = 'energy'
[]
[wall_pressure_x]
type = PCNSFVImplicitMomentumPressureBC
momentum_component = 'x'
boundary = 'left right'
variable = sup_mom_x
[]
[wall_pressure_y]
type = PCNSFVImplicitMomentumPressureBC
momentum_component = 'y'
boundary = 'left right'
variable = sup_mom_y
[]
# Use these to help create more accurate cell centered gradients for cells adjacent to boundaries
[T_bottom]
type = FVDirichletBC
variable = T_fluid
value = ${T}
boundary = 'bottom'
[]
[sup_mom_x_bottom_and_walls]
type = FVDirichletBC
variable = sup_mom_x
value = 0
boundary = 'bottom left right'
[]
[sup_mom_y_walls]
type = FVDirichletBC
variable = sup_mom_y
value = 0
boundary = 'left right'
[]
[sup_mom_y_bottom]
type = FVDirichletBC
variable = sup_mom_y
value = ${sup_mom_y_in}
boundary = 'bottom'
[]
[p_top]
type = FVDirichletBC
variable = pressure
value = ${p_initial}
boundary = 'top'
[]
[]
[Functions]
[ud_in]
type = ParsedVectorFunction
expression_x = '0'
expression_y = '${sup_mom_y_in}'
[]
[eps]
type = ParsedFunction
expression = 'if(y < 2.8, 1,
if(y < 3.2, 1 - .5 / .4 * (y - 2.8),
if(y < 6.8, .5,
if(y < 7.2, .5 - .25 / .4 * (y - 6.8),
if(y < 10.8, .25,
if(y < 11.2, .25 + .25 / .4 * (y - 10.8),
if(y < 14.8, .5,
if(y < 15.2, .5 + .5 / .4 * (y - 14.8),
1))))))))'
[]
[]
[Materials]
[var_mat]
type = PorousMixedVarMaterial
pressure = pressure
T_fluid = T_fluid
superficial_rhou = sup_mom_x
superficial_rhov = sup_mom_y
porosity = porosity
[]
[porosity]
type = GenericFunctionMaterial
prop_names = 'porosity'
prop_values = 'eps'
[]
[ad_generic]
type = ADGenericConstantVectorMaterial
prop_names = 'cl'
prop_values = '${friction_coeff} ${friction_coeff} ${friction_coeff}'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
solve_type = NEWTON
line_search = 'bt'
type = Transient
nl_max_its = 20
[TimeStepper]
type = IterationAdaptiveDT
dt = 5e-5
optimal_iterations = 6
growth_factor = 1.2
[]
num_steps = 10
nl_abs_tol = 1e-8
automatic_scaling = true
compute_scaling_once = false
resid_vs_jac_scaling_param = 0.5
verbose = true
steady_state_detection = true
steady_state_tolerance = 1e-8
normalize_solution_diff_norm_by_dt = false
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
[]
[Outputs]
[out]
type = Exodus
[]
[]
[Debug]
show_var_residual_norms = true
[]
[Postprocessors]
active = ''
[num_nl]
type = NumNonlinearIterations
[]
[total_nl]
type = CumulativeValuePostprocessor
postprocessor = num_nl
[]
[]
(modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/basic-conserved-pcnsfv-kt.i)
[GlobalParams]
fp = fp
limiter = 'central_difference'
two_term_boundary_expansion = true
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 1
xmin = .1
xmax = .6
nx = 2
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[Variables]
[rho]
type = MooseVariableFVReal
[]
[rho_ud]
type = MooseVariableFVReal
[]
[rho_et]
type = MooseVariableFVReal
[]
[]
[ICs]
[pressure]
type = FunctionIC
variable = rho
function = 'exact_rho'
[]
[sup_vel_x]
type = FunctionIC
variable = rho_ud
function = 'exact_rho_ud'
[]
[T_fluid]
type = FunctionIC
variable = rho_et
function = 'exact_rho_et'
[]
[]
[FVKernels]
[mass_advection]
type = PCNSFVKT
variable = rho
eqn = "mass"
[]
[mass_fn]
type = FVBodyForce
variable = rho
function = 'forcing_rho'
[]
[momentum_x_advection]
type = PCNSFVKT
variable = rho_ud
momentum_component = x
eqn = "momentum"
[]
[momentum_fn]
type = FVBodyForce
variable = rho_ud
function = 'forcing_rho_ud'
[]
[fluid_energy_advection]
type = PCNSFVKT
variable = rho_et
eqn = "energy"
[]
[energy_fn]
type = FVBodyForce
variable = rho_et
function = 'forcing_rho_et'
[]
[]
[FVBCs]
[mass_left]
variable = rho
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'mass'
[]
[momentum_left]
variable = rho_ud
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'momentum'
momentum_component = 'x'
[]
[energy_left]
variable = rho_et
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'energy'
[]
[mass_right]
variable = rho
type = PCNSFVStrongBC
boundary = right
eqn = 'mass'
pressure = 'exact_p'
[]
[momentum_right]
variable = rho_ud
type = PCNSFVStrongBC
boundary = right
eqn = 'momentum'
momentum_component = 'x'
pressure = 'exact_p'
[]
[energy_right]
variable = rho_et
type = PCNSFVStrongBC
boundary = right
eqn = 'energy'
pressure = 'exact_p'
[]
# help gradient reconstruction
[rho_right]
type = FVFunctionDirichletBC
variable = rho
function = exact_rho
boundary = 'right'
[]
[rho_ud_left]
type = FVFunctionDirichletBC
variable = rho_ud
function = exact_rho_ud
boundary = 'left'
[]
[rho_et_left]
type = FVFunctionDirichletBC
variable = rho_et
function = exact_rho_et
boundary = 'left'
[]
[]
[Materials]
[var_mat]
type = PorousConservedVarMaterial
rho = rho
superficial_rhou = rho_ud
rho_et = rho_et
porosity = porosity
[]
[porosity]
type = GenericFunctionMaterial
prop_names = 'porosity'
prop_values = 'eps'
[]
[]
[Functions]
[exact_rho]
type = ParsedFunction
expression = '3.48788261470924*cos(x)'
[]
[forcing_rho]
type = ParsedFunction
expression = '-3.45300378856215*sin(1.1*x)'
[]
[exact_rho_ud]
type = ParsedFunction
expression = '3.13909435323832*cos(1.1*x)'
[]
[forcing_rho_ud]
type = ParsedFunction
expression = '-0.9*(10.6975765229419*cos(1.2*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + 0.9*(10.6975765229419*sin(x)*cos(1.2*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 12.8370918275302*sin(1.2*x)/cos(x))*cos(x) + 3.13909435323832*sin(x)*cos(1.1*x)^2/cos(x)^2 - 6.9060075771243*sin(1.1*x)*cos(1.1*x)/cos(x)'
[]
[exact_rho_et]
type = ParsedFunction
expression = '26.7439413073546*cos(1.2*x)'
[]
[forcing_rho_et]
type = ParsedFunction
expression = '0.9*(3.48788261470924*(3.06706896551724*cos(1.2*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.2*x))*sin(x)*cos(1.1*x)/cos(x)^2 - 0.99*(3.48788261470924*(3.06706896551724*cos(1.2*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.2*x))*sin(1.1*x)/cos(x) + 0.9*(-(10.6975765229419*cos(1.2*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + (10.6975765229419*sin(x)*cos(1.2*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 12.8370918275302*sin(1.2*x)/cos(x))*cos(x) - 32.0927295688256*sin(1.2*x))*cos(1.1*x)/cos(x)'
[]
[exact_T]
type = ParsedFunction
expression = '0.0106975765229418*cos(1.2*x)/cos(x) - 0.000697576522941848*cos(1.1*x)^2/cos(x)^2'
[]
[exact_eps_p]
type = ParsedFunction
expression = '3.13909435323832*(3.06706896551724*cos(1.2*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)'
[]
[exact_p]
type = ParsedFunction
expression = '3.48788261470924*(3.06706896551724*cos(1.2*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)'
[]
[exact_sup_vel_x]
type = ParsedFunction
expression = '0.9*cos(1.1*x)/cos(x)'
[]
[exact_superficial_velocity]
type = ParsedVectorFunction
expression_x = '0.9*cos(1.1*x)/cos(x)'
[]
[eps]
type = ParsedFunction
expression = '0.9'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
solve_type = NEWTON
type = Transient
num_steps = 1
dtmin = 1
petsc_options = '-snes_linesearch_monitor'
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
nl_max_its = 50
line_search = bt
[]
[Outputs]
exodus = true
csv = true
[]
[Debug]
show_var_residual_norms = true
[]
[Postprocessors]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2rho]
type = ElementL2Error
variable = rho
function = exact_rho
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2rho_ud]
variable = rho_ud
function = exact_rho_ud
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2rho_et]
variable = rho_et
function = exact_rho_et
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/cns/straight_channel_porosity_step/rotated-2d-bkt-function-porosity.i)
p_initial=1.01e5
T=273.15
# u refers to the superficial velocity
u_in=1
user_limiter='upwind'
friction_coeff=10
[GlobalParams]
fp = fp
two_term_boundary_expansion = true
limiter = ${user_limiter}
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1
nx = 3
ymin = 0
ymax = 18
ny = 90
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[Variables]
[pressure]
type = MooseVariableFVReal
initial_condition = ${p_initial}
[]
[sup_vel_x]
type = MooseVariableFVReal
initial_condition = 1e-15
scaling = 1e-2
[]
[sup_vel_y]
type = MooseVariableFVReal
initial_condition = 1e-15
scaling = 1e-2
[]
[T_fluid]
type = MooseVariableFVReal
initial_condition = ${T}
scaling = 1e-5
[]
[]
[AuxVariables]
[vel_y]
type = MooseVariableFVReal
[]
[sup_mom_y]
type = MooseVariableFVReal
[]
[rho]
type = MooseVariableFVReal
[]
[eps]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[vel_y]
type = ADMaterialRealAux
variable = vel_y
property = vel_y
execute_on = 'timestep_end'
[]
[sup_mom_y]
type = ADMaterialRealAux
variable = sup_mom_y
property = superficial_rhov
execute_on = 'timestep_end'
[]
[rho]
type = ADMaterialRealAux
variable = rho
property = rho
execute_on = 'timestep_end'
[]
[eps]
type = MaterialRealAux
variable = eps
property = porosity
execute_on = 'timestep_end'
[]
[]
[FVKernels]
[mass_time]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rho_dt'
variable = pressure
[]
[mass_advection]
type = PCNSFVKT
variable = pressure
eqn = "mass"
[]
[momentum_time]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rhou_dt'
variable = sup_vel_x
[]
[momentum_advection]
type = PCNSFVKT
variable = sup_vel_x
eqn = "momentum"
momentum_component = 'x'
[]
[eps_grad]
type = PNSFVPGradEpsilon
variable = sup_vel_x
momentum_component = 'x'
epsilon_function = 'eps'
[]
[drag]
type = PCNSFVMomentumFriction
variable = sup_vel_x
momentum_component = 'x'
Darcy_name = 'cl'
momentum_name = superficial_rhou
[]
[momentum_time_y]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rhov_dt'
variable = sup_vel_y
[]
[momentum_advection_y]
type = PCNSFVKT
variable = sup_vel_y
eqn = "momentum"
momentum_component = 'y'
[]
[eps_grad_y]
type = PNSFVPGradEpsilon
variable = sup_vel_y
momentum_component = 'y'
epsilon_function = 'eps'
[]
[drag_y]
type = PCNSFVMomentumFriction
variable = sup_vel_y
momentum_component = 'y'
Darcy_name = 'cl'
momentum_name = superficial_rhov
[]
[energy_time]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rho_et_dt'
variable = T_fluid
[]
[energy_advection]
type = PCNSFVKT
variable = T_fluid
eqn = "energy"
[]
[]
[FVBCs]
[rho_bottom]
type = PCNSFVStrongBC
boundary = 'bottom'
variable = pressure
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'mass'
[]
[rhou_bottom]
type = PCNSFVStrongBC
boundary = 'bottom'
variable = sup_vel_x
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'momentum'
momentum_component = 'x'
[]
[rhov_bottom]
type = PCNSFVStrongBC
boundary = 'bottom'
variable = sup_vel_y
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'momentum'
momentum_component = 'y'
[]
[rho_et_bottom]
type = PCNSFVStrongBC
boundary = 'bottom'
variable = T_fluid
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'energy'
[]
[rho_top]
type = PCNSFVStrongBC
boundary = 'top'
variable = pressure
pressure = ${p_initial}
eqn = 'mass'
[]
[rhou_top]
type = PCNSFVStrongBC
boundary = 'top'
variable = sup_vel_x
pressure = ${p_initial}
eqn = 'momentum'
momentum_component = 'x'
[]
[rhov_top]
type = PCNSFVStrongBC
boundary = 'top'
variable = sup_vel_y
pressure = ${p_initial}
eqn = 'momentum'
momentum_component = 'y'
[]
[rho_et_top]
type = PCNSFVStrongBC
boundary = 'top'
variable = T_fluid
pressure = ${p_initial}
eqn = 'energy'
[]
[wall_pressure_x]
type = PCNSFVImplicitMomentumPressureBC
momentum_component = 'x'
boundary = 'left right'
variable = sup_vel_x
[]
[wall_pressure_y]
type = PCNSFVImplicitMomentumPressureBC
momentum_component = 'y'
boundary = 'left right'
variable = sup_vel_y
[]
# Use these to help create more accurate cell centered gradients for cells adjacent to boundaries
[T_bottom]
type = FVDirichletBC
variable = T_fluid
value = ${T}
boundary = 'bottom'
[]
[sup_vel_x_bottom_and_walls]
type = FVDirichletBC
variable = sup_vel_x
value = 0
boundary = 'bottom left right'
[]
[sup_vel_y_walls]
type = FVDirichletBC
variable = sup_vel_y
value = 0
boundary = 'left right'
[]
[sup_vel_y_bottom]
type = FVDirichletBC
variable = sup_vel_y
value = ${u_in}
boundary = 'bottom'
[]
[p_top]
type = FVDirichletBC
variable = pressure
value = ${p_initial}
boundary = 'top'
[]
[]
[Functions]
[ud_in]
type = ParsedVectorFunction
expression_x = '0'
expression_y = '${u_in}'
[]
[eps]
type = ParsedFunction
expression = 'if(y < 2.8, 1,
if(y < 3.2, 1 - .5 / .4 * (y - 2.8),
if(y < 6.8, .5,
if(y < 7.2, .5 - .25 / .4 * (y - 6.8),
if(y < 10.8, .25,
if(y < 11.2, .25 + .25 / .4 * (y - 10.8),
if(y < 14.8, .5,
if(y < 15.2, .5 + .5 / .4 * (y - 14.8),
1))))))))'
[]
[]
[Materials]
[var_mat]
type = PorousPrimitiveVarMaterial
pressure = pressure
T_fluid = T_fluid
superficial_vel_x = sup_vel_x
superficial_vel_y = sup_vel_y
fp = fp
porosity = porosity
[]
[porosity]
type = GenericFunctionMaterial
prop_names = 'porosity'
prop_values = 'eps'
[]
[ad_generic]
type = ADGenericConstantVectorMaterial
prop_names = 'cl'
prop_values = '${friction_coeff} ${friction_coeff} ${friction_coeff}'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
solve_type = NEWTON
line_search = 'bt'
type = Transient
nl_max_its = 20
[TimeStepper]
type = IterationAdaptiveDT
dt = 5e-5
optimal_iterations = 6
growth_factor = 1.2
[]
num_steps = 10000
end_time = 500
nl_abs_tol = 1e-7
petsc_options_iname = '-pc_type -pc_factor_mat_solver_type'
petsc_options_value = 'lu mumps'
[]
[Outputs]
[out]
type = Exodus
execute_on = 'final'
[]
checkpoint = true
[]
[Debug]
show_var_residual_norms = true
[]
(modules/navier_stokes/test/tests/finite_volume/cns/straight_channel_porosity_step/hllc.i)
p_initial=1.01e5
T=273.15
# u refers to the superficial velocity
u_in=1
[GlobalParams]
fp = fp
two_term_boundary_expansion = true
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 1
xmin = 0
xmax = 18
nx = 180
[]
[to_pt5]
input = cartesian
type = SubdomainBoundingBoxGenerator
bottom_left = '2 0 0'
top_right = '4 1 0'
block_id = 1
[]
[pt5]
input = to_pt5
type = SubdomainBoundingBoxGenerator
bottom_left = '4 0 0'
top_right = '6 1 0'
block_id = 2
[]
[to_pt25]
input = pt5
type = SubdomainBoundingBoxGenerator
bottom_left = '6 0 0'
top_right = '8 1 0'
block_id = 3
[]
[pt25]
input = to_pt25
type = SubdomainBoundingBoxGenerator
bottom_left = '8 0 0'
top_right = '10 1 0'
block_id = 4
[]
[to_pt5_again]
input = pt25
type = SubdomainBoundingBoxGenerator
bottom_left = '10 0 0'
top_right = '12 1 0'
block_id = 5
[]
[pt5_again]
input = to_pt5_again
type = SubdomainBoundingBoxGenerator
bottom_left = '12 0 0'
top_right = '14 1 0'
block_id = 6
[]
[to_one]
input = pt5_again
type = SubdomainBoundingBoxGenerator
bottom_left = '14 0 0'
top_right = '16 1 0'
block_id = 7
[]
[one]
input = to_one
type = SubdomainBoundingBoxGenerator
bottom_left = '16 0 0'
top_right = '18 1 0'
block_id = 8
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[Variables]
[pressure]
type = MooseVariableFVReal
initial_condition = ${p_initial}
[]
[sup_vel_x]
type = MooseVariableFVReal
initial_condition = 1
scaling = 1e-2
[]
[T_fluid]
type = MooseVariableFVReal
initial_condition = ${T}
scaling = 1e-5
[]
[]
[AuxVariables]
[vel_x]
type = MooseVariableFVReal
[]
[sup_mom_x]
type = MooseVariableFVReal
[]
[rho]
type = MooseVariableFVReal
[]
[worst_courant]
type = MooseVariableFVReal
[]
[porosity]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[vel_x]
type = ADMaterialRealAux
variable = vel_x
property = vel_x
execute_on = 'timestep_end'
[]
[sup_mom_x]
type = ADMaterialRealAux
variable = sup_mom_x
property = superficial_rhou
execute_on = 'timestep_end'
[]
[rho]
type = ADMaterialRealAux
variable = rho
property = rho
execute_on = 'timestep_end'
[]
[worst_courant]
type = Courant
variable = worst_courant
u = sup_vel_x
execute_on = 'timestep_end'
[]
[porosity]
type = MaterialRealAux
variable = porosity
property = porosity
execute_on = 'timestep_end'
[]
[]
[FVKernels]
[mass_advection]
type = PCNSFVMassHLLC
variable = pressure
[]
[momentum_advection]
type = PCNSFVMomentumHLLC
variable = sup_vel_x
momentum_component = 'x'
[]
[eps_grad]
type = PNSFVPGradEpsilon
variable = sup_vel_x
momentum_component = 'x'
epsilon_function = 'eps'
[]
[energy_advection]
type = PCNSFVFluidEnergyHLLC
variable = T_fluid
[]
[]
[FVBCs]
[rho_left]
type = PCNSFVStrongBC
boundary = 'left'
variable = pressure
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'mass'
[]
[rhou_left]
type = PCNSFVStrongBC
boundary = 'left'
variable = sup_vel_x
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'momentum'
momentum_component = 'x'
[]
[rho_et_left]
type = PCNSFVStrongBC
boundary = 'left'
variable = T_fluid
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'energy'
[]
[rho_right]
type = PCNSFVStrongBC
boundary = 'right'
variable = pressure
pressure = ${p_initial}
eqn = 'mass'
[]
[rhou_right]
type = PCNSFVStrongBC
boundary = 'right'
variable = sup_vel_x
pressure = ${p_initial}
eqn = 'momentum'
momentum_component = 'x'
[]
[rho_et_right]
type = PCNSFVStrongBC
boundary = 'right'
variable = T_fluid
pressure = ${p_initial}
eqn = 'energy'
[]
# Use these to help create more accurate cell centered gradients for cells adjacent to boundaries
[T_left]
type = FVDirichletBC
variable = T_fluid
value = ${T}
boundary = 'left'
[]
[sup_vel_left]
type = FVDirichletBC
variable = sup_vel_x
value = ${u_in}
boundary = 'left'
[]
[p_right]
type = FVDirichletBC
variable = pressure
value = ${p_initial}
boundary = 'right'
[]
[]
[Functions]
[ud_in]
type = ParsedVectorFunction
expression_x = '${u_in}'
[]
[eps]
type = ParsedFunction
expression = 'if(x < 2, 1,
if(x < 4, 1 - .5 / 2 * (x - 2),
if(x < 6, .5,
if(x < 8, .5 - .25 / 2 * (x - 6),
if(x < 10, .25,
if(x < 12, .25 + .25 / 2 * (x - 10),
if(x < 14, .5,
if(x < 16, .5 + .5 / 2 * (x - 14),
1))))))))'
[]
[]
[Materials]
[var_mat]
type = PorousPrimitiveVarMaterial
pressure = pressure
T_fluid = T_fluid
superficial_vel_x = sup_vel_x
porosity = porosity
[]
[porosity]
type = GenericFunctionMaterial
prop_names = 'porosity'
prop_values = 'eps'
[]
[]
[Executioner]
solve_type = NEWTON
line_search = 'bt'
type = Steady
[]
[Outputs]
[out]
type = Exodus
execute_on = 'final'
[]
checkpoint = true
[]
[Debug]
show_var_residual_norms = true
[]
(modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/varying-eps-hllc.i)
[GlobalParams]
fp = fp
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 1
xmin = .1
xmax = .6
nx = 2
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[Variables]
[pressure]
type = MooseVariableFVReal
[]
[sup_mom_x]
type = MooseVariableFVReal
[]
[T_fluid]
type = MooseVariableFVReal
[]
[]
[ICs]
[pressure]
type = FunctionIC
variable = pressure
function = 'exact_p'
[]
[sup_mom_x]
type = FunctionIC
variable = sup_mom_x
function = 'exact_rho_ud'
[]
[T_fluid]
type = FunctionIC
variable = T_fluid
function = 'exact_T'
[]
[]
[FVKernels]
[mass_advection]
type = PCNSFVMassHLLC
variable = pressure
[]
[mass_fn]
type = FVBodyForce
variable = pressure
function = 'forcing_rho'
[]
[momentum_x_advection]
type = PCNSFVMomentumHLLC
variable = sup_mom_x
momentum_component = x
[]
[eps_grad]
type = PNSFVPGradEpsilon
variable = sup_mom_x
momentum_component = 'x'
epsilon_function = 'eps'
[]
[momentum_fn]
type = FVBodyForce
variable = sup_mom_x
function = 'forcing_rho_ud'
[]
[fluid_energy_advection]
type = PCNSFVFluidEnergyHLLC
variable = T_fluid
[]
[energy_fn]
type = FVBodyForce
variable = T_fluid
function = 'forcing_rho_et'
[]
[]
[FVBCs]
[mass_left]
variable = pressure
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'mass'
[]
[momentum_left]
variable = sup_mom_x
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'momentum'
momentum_component = 'x'
[]
[energy_left]
variable = T_fluid
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'energy'
[]
[mass_right]
variable = pressure
type = PCNSFVStrongBC
boundary = right
eqn = 'mass'
pressure = 'exact_p'
[]
[momentum_right]
variable = sup_mom_x
type = PCNSFVStrongBC
boundary = right
eqn = 'momentum'
momentum_component = 'x'
pressure = 'exact_p'
[]
[energy_right]
variable = T_fluid
type = PCNSFVStrongBC
boundary = right
eqn = 'energy'
pressure = 'exact_p'
[]
[]
[Materials]
[var_mat]
type = PorousMixedVarMaterial
pressure = pressure
superficial_rhou = sup_mom_x
T_fluid = T_fluid
porosity = porosity
[]
[porosity]
type = GenericFunctionMaterial
prop_names = 'porosity'
prop_values = 'eps'
[]
[]
[Functions]
[exact_rho]
type = ParsedFunction
expression = '3.48788261470924*cos(x)'
[]
[forcing_rho]
type = ParsedFunction
expression = '-3.83667087618017*sin(1.1*x)*cos(1.3*x) - 4.53424739912202*sin(1.3*x)*cos(1.1*x)'
[]
[exact_rho_ud]
type = ParsedFunction
expression = '3.48788261470924*cos(1.1*x)*cos(1.3*x)'
[]
[forcing_rho_ud]
type = ParsedFunction
expression = '(-(10.6975765229419*cos(1.5*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + (10.6975765229419*sin(x)*cos(1.5*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 16.0463647844128*sin(1.5*x)/cos(x))*cos(x))*cos(1.3*x) + 3.48788261470924*sin(x)*cos(1.1*x)^2*cos(1.3*x)/cos(x)^2 - 7.67334175236034*sin(1.1*x)*cos(1.1*x)*cos(1.3*x)/cos(x) - 4.53424739912202*sin(1.3*x)*cos(1.1*x)^2/cos(x)'
[]
[exact_rho_et]
type = ParsedFunction
expression = '26.7439413073546*cos(1.5*x)'
[]
[forcing_rho_et]
type = ParsedFunction
expression = '1.0*(3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.5*x))*sin(x)*cos(1.1*x)*cos(1.3*x)/cos(x)^2 - 1.1*(3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.5*x))*sin(1.1*x)*cos(1.3*x)/cos(x) - 1.3*(3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.5*x))*sin(1.3*x)*cos(1.1*x)/cos(x) + 1.0*(-(10.6975765229419*cos(1.5*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + (10.6975765229419*sin(x)*cos(1.5*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 16.0463647844128*sin(1.5*x)/cos(x))*cos(x) - 40.1159119610319*sin(1.5*x))*cos(1.1*x)*cos(1.3*x)/cos(x)'
[]
[exact_T]
type = ParsedFunction
expression = '0.0106975765229418*cos(1.5*x)/cos(x) - 0.000697576522941848*cos(1.1*x)^2/cos(x)^2'
[]
[exact_eps_p]
type = ParsedFunction
expression = '3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)*cos(1.3*x)'
[]
[exact_p]
type = ParsedFunction
expression = '3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)'
[]
[exact_sup_vel_x]
type = ParsedFunction
expression = '1.0*cos(1.1*x)*cos(1.3*x)/cos(x)'
[]
[eps]
type = ParsedFunction
expression = 'cos(1.3*x)'
[]
[exact_superficial_velocity]
type = ParsedVectorFunction
expression_x = '1.0*cos(1.1*x)*cos(1.3*x)/cos(x)'
[]
[]
[Executioner]
solve_type = NEWTON
type = Steady
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
nl_max_its = 50
line_search = bt
nl_rel_tol = 1e-12
nl_abs_tol = 1e-12
[]
[Outputs]
exodus = true
csv = true
[]
[Debug]
show_var_residual_norms = true
[]
[Postprocessors]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2pressure]
type = ElementL2Error
variable = pressure
function = exact_p
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2sup_mom_x]
variable = sup_mom_x
function = exact_rho_ud
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2T_fluid]
variable = T_fluid
function = exact_T
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/cns/straight_channel_porosity_step/implicit-euler-basic-kt-primitive.i)
p_initial=1.01e5
T=273.15
# u refers to the superficial velocity
u_in=1
user_limiter='upwind'
[GlobalParams]
fp = fp
two_term_boundary_expansion = true
limiter = ${user_limiter}
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 1
xmin = 0
xmax = 18
nx = 180
[]
[to_pt5]
input = cartesian
type = SubdomainBoundingBoxGenerator
bottom_left = '2 0 0'
top_right = '4 1 0'
block_id = 1
[]
[pt5]
input = to_pt5
type = SubdomainBoundingBoxGenerator
bottom_left = '4 0 0'
top_right = '6 1 0'
block_id = 2
[]
[to_pt25]
input = pt5
type = SubdomainBoundingBoxGenerator
bottom_left = '6 0 0'
top_right = '8 1 0'
block_id = 3
[]
[pt25]
input = to_pt25
type = SubdomainBoundingBoxGenerator
bottom_left = '8 0 0'
top_right = '10 1 0'
block_id = 4
[]
[to_pt5_again]
input = pt25
type = SubdomainBoundingBoxGenerator
bottom_left = '10 0 0'
top_right = '12 1 0'
block_id = 5
[]
[pt5_again]
input = to_pt5_again
type = SubdomainBoundingBoxGenerator
bottom_left = '12 0 0'
top_right = '14 1 0'
block_id = 6
[]
[to_one]
input = pt5_again
type = SubdomainBoundingBoxGenerator
bottom_left = '14 0 0'
top_right = '16 1 0'
block_id = 7
[]
[one]
input = to_one
type = SubdomainBoundingBoxGenerator
bottom_left = '16 0 0'
top_right = '18 1 0'
block_id = 8
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[Variables]
[pressure]
type = MooseVariableFVReal
initial_condition = ${p_initial}
[]
[sup_vel_x]
type = MooseVariableFVReal
initial_condition = 1e-15
scaling = 1e-2
[]
[T_fluid]
type = MooseVariableFVReal
initial_condition = ${T}
scaling = 1e-5
[]
[]
[AuxVariables]
[vel_x]
type = MooseVariableFVReal
[]
[sup_mom_x]
type = MooseVariableFVReal
[]
[rho]
type = MooseVariableFVReal
[]
[worst_courant]
type = MooseVariableFVReal
[]
[porosity]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[vel_x]
type = ADMaterialRealAux
variable = vel_x
property = vel_x
execute_on = 'timestep_end'
[]
[sup_mom_x]
type = ADMaterialRealAux
variable = sup_mom_x
property = superficial_rhou
execute_on = 'timestep_end'
[]
[rho]
type = ADMaterialRealAux
variable = rho
property = rho
execute_on = 'timestep_end'
[]
[worst_courant]
type = Courant
variable = worst_courant
u = sup_vel_x
execute_on = 'timestep_end'
[]
[porosity]
type = MaterialRealAux
variable = porosity
property = porosity
execute_on = 'timestep_end'
[]
[]
[FVKernels]
[mass_time]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rho_dt'
variable = pressure
[]
[mass_advection]
type = PCNSFVKT
variable = pressure
eqn = "mass"
[]
[momentum_time]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rhou_dt'
variable = sup_vel_x
[]
[momentum_advection]
type = PCNSFVKT
variable = sup_vel_x
eqn = "momentum"
momentum_component = 'x'
[]
[eps_grad]
type = PNSFVPGradEpsilon
variable = sup_vel_x
momentum_component = 'x'
epsilon_function = 'eps'
[]
[energy_time]
type = FVMatPropTimeKernel
mat_prop_time_derivative = 'dsuperficial_rho_et_dt'
variable = T_fluid
[]
[energy_advection]
type = PCNSFVKT
variable = T_fluid
eqn = "energy"
[]
[]
[FVBCs]
[rho_left]
type = PCNSFVStrongBC
boundary = 'left'
variable = pressure
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'mass'
[]
[rhou_left]
type = PCNSFVStrongBC
boundary = 'left'
variable = sup_vel_x
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'momentum'
momentum_component = 'x'
[]
[rho_et_left]
type = PCNSFVStrongBC
boundary = 'left'
variable = T_fluid
superficial_velocity = 'ud_in'
T_fluid = ${T}
eqn = 'energy'
[]
[rho_right]
type = PCNSFVStrongBC
boundary = 'right'
variable = pressure
pressure = ${p_initial}
eqn = 'mass'
[]
[rhou_right]
type = PCNSFVStrongBC
boundary = 'right'
variable = sup_vel_x
pressure = ${p_initial}
eqn = 'momentum'
momentum_component = 'x'
[]
[rho_et_right]
type = PCNSFVStrongBC
boundary = 'right'
variable = T_fluid
pressure = ${p_initial}
eqn = 'energy'
[]
# Use these to help create more accurate cell centered gradients for cells adjacent to boundaries
[T_left]
type = FVDirichletBC
variable = T_fluid
value = ${T}
boundary = 'left'
[]
[sup_vel_left]
type = FVDirichletBC
variable = sup_vel_x
value = ${u_in}
boundary = 'left'
[]
[p_right]
type = FVDirichletBC
variable = pressure
value = ${p_initial}
boundary = 'right'
[]
[]
[Functions]
[ud_in]
type = ParsedVectorFunction
expression_x = '${u_in}'
[]
[eps]
type = ParsedFunction
expression = 'if(x < 2, 1,
if(x < 4, 1 - .5 / 2 * (x - 2),
if(x < 6, .5,
if(x < 8, .5 - .25 / 2 * (x - 6),
if(x < 10, .25,
if(x < 12, .25 + .25 / 2 * (x - 10),
if(x < 14, .5,
if(x < 16, .5 + .5 / 2 * (x - 14),
1))))))))'
[]
[]
[Materials]
[var_mat]
type = PorousPrimitiveVarMaterial
pressure = pressure
T_fluid = T_fluid
superficial_vel_x = sup_vel_x
fp = fp
porosity = porosity
[]
[porosity]
type = GenericFunctionMaterial
prop_names = 'porosity'
prop_values = 'eps'
[]
[]
[Executioner]
solve_type = NEWTON
line_search = 'bt'
type = Transient
nl_max_its = 20
[TimeStepper]
type = IterationAdaptiveDT
dt = 5e-5
optimal_iterations = 6
growth_factor = 1.2
[]
num_steps = 10000
end_time = 500
nl_abs_tol = 1e-8
[]
[Outputs]
[out]
type = Exodus
execute_on = 'final'
[]
checkpoint = true
[]
[Debug]
show_var_residual_norms = true
[]
(modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/varying-eps-basic-kt-mixed.i)
[GlobalParams]
fp = fp
limiter = 'central_difference'
two_term_boundary_expansion = true
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 1
xmin = .1
xmax = .6
nx = 2
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[Variables]
[pressure]
type = MooseVariableFVReal
[]
[sup_mom_x]
type = MooseVariableFVReal
[]
[T_fluid]
type = MooseVariableFVReal
[]
[]
[ICs]
[pressure]
type = FunctionIC
variable = pressure
function = 'exact_p'
[]
[sup_mom_x]
type = FunctionIC
variable = sup_mom_x
function = 'exact_rho_ud'
[]
[T_fluid]
type = FunctionIC
variable = T_fluid
function = 'exact_T'
[]
[]
[FVKernels]
[mass_advection]
type = PCNSFVKT
variable = pressure
eqn = "mass"
[]
[mass_fn]
type = FVBodyForce
variable = pressure
function = 'forcing_rho'
[]
[momentum_x_advection]
type = PCNSFVKT
variable = sup_mom_x
momentum_component = x
eqn = "momentum"
[]
[eps_grad]
type = PNSFVPGradEpsilon
variable = sup_mom_x
momentum_component = 'x'
epsilon_function = 'eps'
[]
[momentum_fn]
type = FVBodyForce
variable = sup_mom_x
function = 'forcing_rho_ud'
[]
[fluid_energy_advection]
type = PCNSFVKT
variable = T_fluid
eqn = "energy"
[]
[energy_fn]
type = FVBodyForce
variable = T_fluid
function = 'forcing_rho_et'
[]
[]
[FVBCs]
[mass_left]
variable = pressure
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'mass'
[]
[momentum_left]
variable = sup_mom_x
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'momentum'
momentum_component = 'x'
[]
[energy_left]
variable = T_fluid
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'energy'
[]
[mass_right]
variable = pressure
type = PCNSFVStrongBC
boundary = right
eqn = 'mass'
pressure = 'exact_p'
[]
[momentum_right]
variable = sup_mom_x
type = PCNSFVStrongBC
boundary = right
eqn = 'momentum'
momentum_component = 'x'
pressure = 'exact_p'
[]
[energy_right]
variable = T_fluid
type = PCNSFVStrongBC
boundary = right
eqn = 'energy'
pressure = 'exact_p'
[]
# help gradient reconstruction
[pressure_right]
type = FVFunctionDirichletBC
variable = pressure
function = exact_p
boundary = 'right'
[]
[sup_mom_x_left]
type = FVFunctionDirichletBC
variable = sup_mom_x
function = exact_rho_ud
boundary = 'left'
[]
[T_fluid_left]
type = FVFunctionDirichletBC
variable = T_fluid
function = exact_T
boundary = 'left'
[]
[]
[Materials]
[var_mat]
type = PorousMixedVarMaterial
pressure = pressure
superficial_rhou = sup_mom_x
T_fluid = T_fluid
porosity = porosity
[]
[porosity]
type = GenericFunctionMaterial
prop_names = 'porosity'
prop_values = 'eps'
[]
[]
[Functions]
[exact_rho]
type = ParsedFunction
expression = '3.48788261470924*cos(x)'
[]
[forcing_rho]
type = ParsedFunction
expression = '-3.83667087618017*sin(1.1*x)*cos(1.3*x) - 4.53424739912202*sin(1.3*x)*cos(1.1*x)'
[]
[exact_rho_ud]
type = ParsedFunction
expression = '3.48788261470924*cos(1.1*x)*cos(1.3*x)'
[]
[forcing_rho_ud]
type = ParsedFunction
expression = '(-(10.6975765229419*cos(1.5*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + (10.6975765229419*sin(x)*cos(1.5*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 16.0463647844128*sin(1.5*x)/cos(x))*cos(x))*cos(1.3*x) + 3.48788261470924*sin(x)*cos(1.1*x)^2*cos(1.3*x)/cos(x)^2 - 7.67334175236034*sin(1.1*x)*cos(1.1*x)*cos(1.3*x)/cos(x) - 4.53424739912202*sin(1.3*x)*cos(1.1*x)^2/cos(x)'
[]
[exact_rho_et]
type = ParsedFunction
expression = '26.7439413073546*cos(1.5*x)'
[]
[forcing_rho_et]
type = ParsedFunction
expression = '1.0*(3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.5*x))*sin(x)*cos(1.1*x)*cos(1.3*x)/cos(x)^2 - 1.1*(3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.5*x))*sin(1.1*x)*cos(1.3*x)/cos(x) - 1.3*(3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.5*x))*sin(1.3*x)*cos(1.1*x)/cos(x) + 1.0*(-(10.6975765229419*cos(1.5*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + (10.6975765229419*sin(x)*cos(1.5*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 16.0463647844128*sin(1.5*x)/cos(x))*cos(x) - 40.1159119610319*sin(1.5*x))*cos(1.1*x)*cos(1.3*x)/cos(x)'
[]
[exact_T]
type = ParsedFunction
expression = '0.0106975765229418*cos(1.5*x)/cos(x) - 0.000697576522941848*cos(1.1*x)^2/cos(x)^2'
[]
[exact_eps_p]
type = ParsedFunction
expression = '3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)*cos(1.3*x)'
[]
[exact_p]
type = ParsedFunction
expression = '3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)'
[]
[exact_sup_vel_x]
type = ParsedFunction
expression = '1.0*cos(1.1*x)*cos(1.3*x)/cos(x)'
[]
[eps]
type = ParsedFunction
expression = 'cos(1.3*x)'
[]
[exact_superficial_velocity]
type = ParsedVectorFunction
expression_x = '1.0*cos(1.1*x)*cos(1.3*x)/cos(x)'
[]
[]
[Executioner]
solve_type = NEWTON
type = Transient
num_steps = 1
dtmin = 1
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
nl_max_its = 50
line_search = bt
nl_rel_tol = 1e-12
nl_abs_tol = 1e-12
[]
[Outputs]
exodus = true
csv = true
[]
[Debug]
show_var_residual_norms = true
[]
[Postprocessors]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2pressure]
type = ElementL2Error
variable = pressure
function = exact_p
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2sup_mom_x]
variable = sup_mom_x
function = exact_rho_ud
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2T_fluid]
variable = T_fluid
function = exact_T
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/basic-primitive-pcnsfv-kt.i)
[GlobalParams]
fp = fp
limiter = 'central_difference'
two_term_boundary_expansion = true
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 1
xmin = .1
xmax = .6
nx = 2
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[Variables]
[pressure]
type = MooseVariableFVReal
[]
[sup_vel_x]
type = MooseVariableFVReal
[]
[T_fluid]
type = MooseVariableFVReal
[]
[]
[ICs]
[pressure]
type = FunctionIC
variable = pressure
function = 'exact_p'
[]
[sup_vel_x]
type = FunctionIC
variable = sup_vel_x
function = 'exact_sup_vel_x'
[]
[T_fluid]
type = FunctionIC
variable = T_fluid
function = 'exact_T'
[]
[]
[FVKernels]
[mass_advection]
type = PCNSFVKT
variable = pressure
eqn = "mass"
[]
[mass_fn]
type = FVBodyForce
variable = pressure
function = 'forcing_rho'
[]
[momentum_x_advection]
type = PCNSFVKT
variable = sup_vel_x
momentum_component = x
eqn = "momentum"
[]
[momentum_fn]
type = FVBodyForce
variable = sup_vel_x
function = 'forcing_rho_ud'
[]
[fluid_energy_advection]
type = PCNSFVKT
variable = T_fluid
eqn = "energy"
[]
[energy_fn]
type = FVBodyForce
variable = T_fluid
function = 'forcing_rho_et'
[]
[]
[FVBCs]
[mass_left]
variable = pressure
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'mass'
[]
[momentum_left]
variable = sup_vel_x
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'momentum'
momentum_component = 'x'
[]
[energy_left]
variable = T_fluid
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'energy'
[]
[mass_right]
variable = pressure
type = PCNSFVStrongBC
boundary = right
eqn = 'mass'
pressure = 'exact_p'
[]
[momentum_right]
variable = sup_vel_x
type = PCNSFVStrongBC
boundary = right
eqn = 'momentum'
momentum_component = 'x'
pressure = 'exact_p'
[]
[energy_right]
variable = T_fluid
type = PCNSFVStrongBC
boundary = right
eqn = 'energy'
pressure = 'exact_p'
[]
# help gradient reconstruction
[pressure_right]
type = FVFunctionDirichletBC
variable = pressure
function = exact_p
boundary = 'right'
[]
[sup_vel_x_left]
type = FVFunctionDirichletBC
variable = sup_vel_x
function = exact_sup_vel_x
boundary = 'left'
[]
[T_fluid_left]
type = FVFunctionDirichletBC
variable = T_fluid
function = exact_T
boundary = 'left'
[]
[]
[Materials]
[var_mat]
type = PorousPrimitiveVarMaterial
pressure = pressure
superficial_vel_x = sup_vel_x
T_fluid = T_fluid
porosity = porosity
[]
[porosity]
type = GenericFunctionMaterial
prop_names = 'porosity'
prop_values = 'eps'
[]
[]
[Functions]
[exact_rho]
type = ParsedFunction
expression = '3.48788261470924*cos(x)'
[]
[forcing_rho]
type = ParsedFunction
expression = '-3.45300378856215*sin(1.1*x)'
[]
[exact_rho_ud]
type = ParsedFunction
expression = '3.13909435323832*cos(1.1*x)'
[]
[forcing_rho_ud]
type = ParsedFunction
expression = '-0.9*(10.6975765229419*cos(1.2*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + 0.9*(10.6975765229419*sin(x)*cos(1.2*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 12.8370918275302*sin(1.2*x)/cos(x))*cos(x) + 3.13909435323832*sin(x)*cos(1.1*x)^2/cos(x)^2 - 6.9060075771243*sin(1.1*x)*cos(1.1*x)/cos(x)'
[]
[exact_rho_et]
type = ParsedFunction
expression = '26.7439413073546*cos(1.2*x)'
[]
[forcing_rho_et]
type = ParsedFunction
expression = '0.9*(3.48788261470924*(3.06706896551724*cos(1.2*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.2*x))*sin(x)*cos(1.1*x)/cos(x)^2 - 0.99*(3.48788261470924*(3.06706896551724*cos(1.2*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.2*x))*sin(1.1*x)/cos(x) + 0.9*(-(10.6975765229419*cos(1.2*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + (10.6975765229419*sin(x)*cos(1.2*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 12.8370918275302*sin(1.2*x)/cos(x))*cos(x) - 32.0927295688256*sin(1.2*x))*cos(1.1*x)/cos(x)'
[]
[exact_T]
type = ParsedFunction
expression = '0.0106975765229418*cos(1.2*x)/cos(x) - 0.000697576522941848*cos(1.1*x)^2/cos(x)^2'
[]
[exact_eps_p]
type = ParsedFunction
expression = '3.13909435323832*(3.06706896551724*cos(1.2*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)'
[]
[exact_p]
type = ParsedFunction
expression = '3.48788261470924*(3.06706896551724*cos(1.2*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)'
[]
[exact_sup_vel_x]
type = ParsedFunction
expression = '0.9*cos(1.1*x)/cos(x)'
[]
[exact_superficial_velocity]
type = ParsedVectorFunction
expression_x = '0.9*cos(1.1*x)/cos(x)'
[]
[eps]
type = ParsedFunction
expression = '0.9'
[]
[]
[Executioner]
solve_type = NEWTON
type = Transient
num_steps = 1
dtmin = 1
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
nl_max_its = 50
line_search = bt
nl_rel_tol = 1e-12
nl_abs_tol = 1e-12
[]
[Outputs]
exodus = true
csv = true
[]
[Debug]
show_var_residual_norms = true
[]
[Postprocessors]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2pressure]
type = ElementL2Error
variable = pressure
function = exact_p
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2sup_vel_x]
variable = sup_vel_x
function = exact_sup_vel_x
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2T_fluid]
variable = T_fluid
function = exact_T
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]