- variableThe name of the variable that this residual object operates on
C++ Type:NonlinearVariableName
Unit:(no unit assumed)
Controllable:No
Description:The name of the variable that this residual object operates on
FVTimeKernel
Residual contribution from time derivative of a variable for the finite volume method.
The time derivative is automatically computed based on the time integration scheme selected. This class should be used in finite volume simulations which leverage the quadrature-point pre-initialized paradigm, which includes fully compressible Navier-Stokes simulations.
When creating a new time derivative kernel, developers should consider inheriting this class as its provides the matrix/vector time tags. If not, those should be added in the validParams()
routine of the new class.
Example input syntax
In this example, the Burger's equation is solved in one dimension. The time derivative term of the equation is added to the numerical system using a FVTimeKernel
.
[FVKernels]
[./burgers]
type = FVBurgers1D
variable = v
[../]
[./time]
type = FVTimeKernel
variable = v
[../]
[]
(test/tests/fvkernels/fv_burgers/fv_burgers.i)Input Parameters
- 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
- 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.
- 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.
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_tagssystem timeThe tag for the matrices this Kernel should fill
Default:system time
C++ Type:MultiMooseEnum
Unit:(no unit assumed)
Options:nontime, system, time
Controllable:No
Description:The tag for the matrices this Kernel should fill
- vector_tagstimeThe tag for the vectors this Kernel should fill
Default:time
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
- seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Unit:(no unit assumed)
Controllable:No
Description:The seed for the master random number generator
- use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Default:False
C++ Type:bool
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
- ghost_layers1The number of layers of elements to ghost.
Default:1
C++ Type:unsigned short
Unit:(no unit assumed)
Controllable:No
Description:The number of layers of elements to ghost.
- use_point_neighborsFalseWhether to use point neighbors, which introduces additional ghosting to that used for simple face neighbors.
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Whether to use point neighbors, which introduces additional ghosting to that used for simple face neighbors.
Parallel Ghosting Parameters
Input Files
- (test/tests/mesh/preparedness/test.i)
- (test/tests/auxkernels/time_derivative_aux/test_fv.i)
- (test/tests/fvkernels/fv_burgers/fv_burgers.i)
- (test/tests/fvkernels/fv_dotdot/fv_dotdot.i)
- (test/tests/fvkernels/constraints/bounded_value.i)
- (test/tests/outputs/debug/show_execution_fv_flux_objects.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/scalar_advection/mass-frac-advection.i)
- (test/tests/indicators/gradient_jump_indicator/gradient_jump_indicator_fv_test.i)
- (test/tests/fvkernels/fv_simple_diffusion/grad-adaptive.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/benchmark_shock_tube_1D/hllc_sod_shocktube.i)
- (test/tests/postprocessors/side_advection_flux_integral/side_advection_flux_integral_fv.i)
- (test/tests/postprocessors/pseudotimestep/fv_burgers_pseudo.i)
- (test/tests/auxkernels/advection_flux/advection_flux_fv.i)
- (test/tests/fvkernels/fv_constant_scalar_advection/2D_constant_scalar_advection.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/implicit_bcs/hllc_sod_shocktube.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/stagnation_inlet/supersonic_nozzle_hllc.i)
- (modules/porous_flow/test/tests/heat_conduction/two_phase_fv.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/shock_tube_2D_cavity/hllc_sod_shocktube_2D.i)
- (test/tests/tag/tag-fv.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/pressure_outlet/subsonic_nozzle_fixed_inflow_hllc.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/symmetry_test/2D_symmetry.i)
Child Objects
- (modules/porous_flow/include/fvkernels/FVPorousFlowMassTimeDerivative.h)
- (modules/navier_stokes/include/fvkernels/FVPorosityTimeDerivative.h)
- (modules/navier_stokes/include/fvkernels/PCNSFVDensityTimeDerivative.h)
- (modules/porous_flow/include/fvkernels/FVPorousFlowEnergyTimeDerivative.h)
- (modules/navier_stokes/include/fvkernels/PCNSFVEnergyTimeDerivative.h)
- (modules/navier_stokes/include/fvkernels/FVMatPropTimeKernel.h)
- (modules/heat_transfer/include/fvkernels/FVHeatConductionTimeDerivative.h)
(test/tests/fvkernels/fv_burgers/fv_burgers.i)
[Mesh]
[./gen_mesh]
type = GeneratedMeshGenerator
dim = 1
xmin = 0
xmax = 10
nx = 50
[../]
[]
[Variables]
[./v]
family = MONOMIAL
order = CONSTANT
fv = true
[../]
[]
[ICs]
[./v_ic]
type = FunctionIC
variable = v
function = 'if (x > 2 & x < 3, 0.5, 0)'
[../]
[]
[FVKernels]
[./burgers]
type = FVBurgers1D
variable = v
[../]
[./time]
type = FVTimeKernel
variable = v
[../]
[]
[FVBCs]
[./fv_burgers_outflow]
type = FVBurgersOutflowBC
variable = v
boundary = 'left right'
[../]
[]
[Executioner]
type = Transient
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
petsc_options = '-snes_converged_reason'
nl_abs_tol = 1e-7
nl_rel_tol = 1e-8
num_steps = 5
dt = 0.05
[]
[Outputs]
exodus = true
[]
(test/tests/mesh/preparedness/test.i)
[GlobalParams]
prevent_boundary_ids_overlap = false
[]
[Mesh]
[region_2_gen]
type = CartesianMeshGenerator
dim = 2
dx = '0.065 0.13 0.305 0.17 0.196'
ix = ' 2 2 2 2 2'
dy = '0.85438 '
iy = '6'
subdomain_id = '68 68 68 68 68'
[]
[region_2_move]
type = TransformGenerator
transform = TRANSLATE
vector_value = '1.2 1.551 0'
input = region_2_gen
[]
[region_3_gen]
type = CartesianMeshGenerator
dim = 2
dx = '0.24 0.24 0.24 0.24 0.24'
ix = ' 2 2 2 2 2'
dy = '0.744166666666666 0.744166666666667 0.744166666666667'
iy = ' 2 2 2'
subdomain_id = '56 57 58 59 60
51 52 53 54 55
46 47 48 49 50'
[]
[region_3_move]
type = TransformGenerator
transform = TRANSLATE
vector_value = '0 2.40538 0'
input = region_3_gen
[]
[region_1_gen]
type = GeneratedMeshGenerator
dim = 2
nx = 10
ny = 6
xmin = 0
xmax = 0.26
ymin = 1.551
ymax = 1.851
subdomain_ids = '62 62 62 62 62 62 62 62 62 62
62 62 62 62 62 62 62 62 62 62
62 62 62 62 62 62 62 62 62 62
62 62 62 62 62 62 62 62 62 62
62 62 62 62 62 62 62 62 62 62
62 62 62 62 62 62 62 62 62 62'
[]
[region_1_extend_1]
type = FillBetweenSidesetsGenerator
input_mesh_1 = 'region_3_move'
input_mesh_2 = 'region_1_gen'
boundary_1 = '0'
boundary_2 = '2'
num_layers = 6
block_id= 61
use_quad_elements = true
keep_inputs = true
begin_side_boundary_id = '3'
end_side_boundary_id = '1'
[]
[region_1_extend_2]
type = FillBetweenSidesetsGenerator
input_mesh_1 = 'region_2_move'
input_mesh_2 = 'region_1_gen'
boundary_1 = 3
boundary_2 = 1
num_layers = 6
block_id= 69
use_quad_elements = true
keep_inputs = false
begin_side_boundary_id = '0'
end_side_boundary_id = '3'
input_boundary_1_id = '1'
input_boundary_2_id = '3'
[]
[region_2_2_gen]
type = CartesianMeshGenerator
dim = 2
dx = '0.065 0.13 0.305 0.17 0.196'
ix = ' 2 2 2 2 2'
dy = '0.85438 '
iy = '6'
subdomain_id = '68 68 68 68 68'
[]
[region_2_2_move]
type = TransformGenerator
transform = TRANSLATE
vector_value = '1.2 1.551 0'
input = region_2_2_gen
[]
[region_6_gen]
type = CartesianMeshGenerator
dim = 2
dx = '0.26 0.94 0.065 0.13 0.305 0.17 0.196'
ix = '10 6 2 2 2 2 2'
dy = '0.584 0.967'
iy = ' 4 6'
subdomain_id = '62 72 72 72 72 72 72
62 70 71 71 71 71 71'
[]
[stitch_1_2_6]
type = StitchedMeshGenerator
inputs = 'region_1_extend_1 region_1_extend_2 region_2_2_move region_6_gen'
stitch_boundaries_pairs = '1 3;
1 3;
0 2'
merge_boundaries_with_same_name = false
[]
[rename_boundary_stitch_1_2_6]
type = RenameBoundaryGenerator
input = stitch_1_2_6
old_boundary = '1'
new_boundary = '2'
[]
[region_4_gen]
type = CartesianMeshGenerator
dim = 2
dx = '0.065 0.13'
ix = ' 2 2 '
dy = '0.744166666666666 0.744166666666667 0.744166666666667'
iy = ' 2 2 2'
subdomain_id = '78 92
78 91
78 90'
[]
[region_4_move]
type = TransformGenerator
transform = TRANSLATE
vector_value = '1.2 2.40538 0'
input = region_4_gen
[]
[region_5_gen]
type = CartesianMeshGenerator
dim = 2
dx = '0.17 0.196'
ix = '2 2'
dy = '0.39 1.8425'
iy = '2 4'
subdomain_id = '100 104
100 104'
[]
[region_5_move]
type = TransformGenerator
transform = TRANSLATE
vector_value = '1.7 2.40538 0'
input = region_5_gen
[]
[region_5_extend]
type = FillBetweenSidesetsGenerator
input_mesh_1 = 'region_4_move'
input_mesh_2 = 'region_5_move'
boundary_1 = 1
boundary_2 = 3
num_layers = 2
block_id= 96
use_quad_elements = true
keep_inputs = true
begin_side_boundary_id = '0'
end_side_boundary_id = '2'
[]
[rename_boundary_region_5]
type = RenameBoundaryGenerator
input = region_5_extend
old_boundary = '0'
new_boundary = '3'
[]
[stitch_1_2_6_5]
type = StitchedMeshGenerator
inputs = 'rename_boundary_stitch_1_2_6 rename_boundary_region_5'
stitch_boundaries_pairs = '2 3;'
merge_boundaries_with_same_name = false
[]
[region_7_gen]
type = CartesianMeshGenerator
dim = 2
dx = '0.24 0.24 0.24 0.24 0.24 0.065 0.13 0.305 0.17 0.196'
ix = ' 2 2 2 2 2 2 2 2 2 2'
dy = '0.744166666666667 0.744166666666667 0.744166666666667 0.744166666666667
0.744166666666667 0.744166666666667 0.744166666666666 0.744166666666666
0.744166666666666 0.458 0.86002'
iy = '2 2 2 2 2 2 2 2 2 2 4'
subdomain_id = '41 42 43 44 45 77 89 95 99 103
36 37 38 39 40 77 88 95 99 103
31 32 33 34 35 77 87 95 99 103
26 27 28 29 30 76 86 94 98 102
21 22 23 24 25 76 85 94 98 102
16 17 18 19 20 76 84 94 98 102
11 12 13 14 15 75 83 93 97 101
6 7 8 9 10 75 82 93 97 101
1 2 3 4 5 75 81 93 97 101
67 67 67 67 67 74 80 65 65 66
63 63 63 63 63 73 79 64 64 64'
[]
[region_7_move]
type = TransformGenerator
transform = TRANSLATE
vector_value = '0.0 4.63788 0'
input = region_7_gen
[]
[stitch]
type = StitchedMeshGenerator
inputs = 'stitch_1_2_6_5 region_7_move'
stitch_boundaries_pairs = '2 0'
merge_boundaries_with_same_name = false
[]
[rename_boundary_1]
type = BoundaryDeletionGenerator
input = stitch
boundary_names = '0 1 2 3'
[]
[rename_boundary_2]
type = SideSetsFromPointsGenerator
input = rename_boundary_1
new_boundary = '2 4 1 3'
points = '1.2 0. 0.
2.066 1.551 0.
1.2 12.6534 0.
0. 1.551 0.'
[]
[rename_boundary_3]
type = RenameBoundaryGenerator
input = rename_boundary_2
new_boundary = 'rbottom ssright rtop ssleft'
old_boundary = '2 4 1 3'
[]
[rename_blocks]
type = RenameBlockGenerator
input = rename_boundary_3
old_block = '1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20
21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40
41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60
61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80
81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100
101 102 103 104'
new_block = 'pbedfuel001 pbedfuel002 pbedfuel003 pbedfuel004 pbedfuel005
pbedfuel006 pbedfuel007 pbedfuel008 pbedfuel009 pbedfuel010
pbedfuel011 pbedfuel012 pbedfuel013 pbedfuel014 pbedfuel015
pbedfuel016 pbedfuel017 pbedfuel018 pbedfuel019 pbedfuel020
pbedfuel021 pbedfuel022 pbedfuel023 pbedfuel024 pbedfuel025
pbedfuel026 pbedfuel027 pbedfuel028 pbedfuel029 pbedfuel030
pbedfuel031 pbedfuel032 pbedfuel033 pbedfuel034 pbedfuel035
pbedfuel036 pbedfuel037 pbedfuel038 pbedfuel039 pbedfuel040
pbedfuel041 pbedfuel042 pbedfuel043 pbedfuel044 pbedfuel045
pbedfuel046 pbedfuel047 pbedfuel048 pbedfuel049 pbedfuel050
pbedfuel051 pbedfuel052 pbedfuel053 pbedfuel054 pbedfuel055
pbedfuel056 pbedfuel057 pbedfuel058 pbedfuel059 pbedfuel060
consfuel061 dischfuel062 upref063 upref064 upref065 upref066
upcvt067 lwref068 outch069 lwrpln070 htleg071 lwref072 buffr073
buffr074 buffr075 buffr076 buffr077 buffr078 crds079 crds080
crds081 crds082 crds083 crds084 crds085 crds086 crds087 crds088
crds089 crds090 crds091 crds092 radrf093 radrf094 radrf095 radrf096
risr097 risr098 risr099 risr100 radrf101 radrf102 radrf103 radrf104'
[]
[]
[Variables]
[T_solid]
type = MooseVariableFVReal
initial_condition = 100
[]
[]
[FVKernels]
[energy_storage]
type = FVTimeKernel
variable = T_solid
[]
[solid_energy_diffusion_core]
type = FVAnisotropicDiffusion
variable = T_solid
coeff = 'effective_thermal_conductivity'
[]
[]
[FVBCs]
[side_set_bc1]
type = FVDirichletBC
variable = T_solid
value = '300'
boundary = 'rtop'
[]
[side_set_bc2]
type = FVDirichletBC
variable = T_solid
value = '600'
boundary = 'rbottom'
[]
[]
[Materials]
[all_channels_porosity]
type = ADGenericFunctorMaterial
prop_names = 'porosity'
prop_values = 0.5
[]
[solid_blocks_full_density_graphite]
type = ADGenericFunctorMaterial
prop_names = 'rho_s cp_s k_s '
prop_values = '1.0 2.0 3.0'
[]
[effective_solid_thermal_conductivity_solid_only]
type = ADGenericVectorFunctorMaterial
prop_names = 'effective_thermal_conductivity'
prop_values = 'k_s k_s k_s'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -ksp_gmres_restart -sub_pc_factor_shift_type'
petsc_options_value = 'lu 100 NONZERO'
# Tolerances.
nl_abs_tol = 1e-8
nl_rel_tol = 1e-9
line_search = none
nl_max_its = 15
[TimeStepper]
type = IterationAdaptiveDT
dt = 0.05
cutback_factor = 0.5
growth_factor = 2.00
optimal_iterations = 6
[]
# Steady state detection.
steady_state_detection = true
steady_state_tolerance = 1e-13
abort_on_solve_fail = true
num_steps = 1
[]
[Outputs]
exodus = true
print_linear_residuals = false
print_linear_converged_reason = false
print_nonlinear_converged_reason = false
[]
(test/tests/auxkernels/time_derivative_aux/test_fv.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 6
ny = 6
[]
[Variables]
[u]
type = MooseVariableFVReal
initial_condition = 2
[]
[]
[FVKernels]
[time]
type = FVTimeKernel
variable = u
[]
[reaction]
type = FVReaction
variable = u
rate = 2.0
[]
[diffusion]
type = FVDiffusion
variable = u
coeff = 0.1
[]
[]
[FVBCs]
[left]
type = FVNeumannBC
variable = u
value = 5
boundary = 'left'
[]
[]
[AuxVariables]
inactive = 'variable_derivative'
[variable_derivative]
family = MONOMIAL
order = CONSTANT
[]
[variable_derivative_fv]
family = MONOMIAL
order = CONSTANT
fv = true
[]
[]
[AuxKernels]
# Time derivative of a FV variable using the functor system
[function_derivative_element]
type = TimeDerivativeAux
variable = variable_derivative_fv
functor = 'u'
factor = 2
[]
# this places the derivative of a FV variable in a FE one
# let's output a warning
inactive = 'function_derivative_element_fv_fe'
[function_derivative_element_fv_fe]
type = TimeDerivativeAux
variable = variable_derivative
functor = 'u'
factor = 2
[]
[]
[Executioner]
type = Transient
dt = 0.1
num_steps = 2
nl_abs_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(test/tests/fvkernels/fv_burgers/fv_burgers.i)
[Mesh]
[./gen_mesh]
type = GeneratedMeshGenerator
dim = 1
xmin = 0
xmax = 10
nx = 50
[../]
[]
[Variables]
[./v]
family = MONOMIAL
order = CONSTANT
fv = true
[../]
[]
[ICs]
[./v_ic]
type = FunctionIC
variable = v
function = 'if (x > 2 & x < 3, 0.5, 0)'
[../]
[]
[FVKernels]
[./burgers]
type = FVBurgers1D
variable = v
[../]
[./time]
type = FVTimeKernel
variable = v
[../]
[]
[FVBCs]
[./fv_burgers_outflow]
type = FVBurgersOutflowBC
variable = v
boundary = 'left right'
[../]
[]
[Executioner]
type = Transient
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
petsc_options = '-snes_converged_reason'
nl_abs_tol = 1e-7
nl_rel_tol = 1e-8
num_steps = 5
dt = 0.05
[]
[Outputs]
exodus = true
[]
(test/tests/fvkernels/fv_dotdot/fv_dotdot.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[v]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 7
[]
[]
[Kernels]
[]
[FVKernels]
[./time]
type = FVTimeKernel
variable = v
[../]
[diff]
type = FVDiffusion
variable = v
coeff = coeff
[]
[]
[FVBCs]
[left]
type = FVDirichletBC
variable = v
boundary = left
value = 7
[]
[right]
type = FVDirichletBC
variable = v
boundary = right
value = 42
[]
[]
[Materials]
[diff]
type = ADGenericFunctorMaterial
prop_names = 'coeff'
prop_values = '.2'
[]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
scheme = newmark-beta
num_steps = 20
dt = 0.1
[]
[Postprocessors]
[vdotdot]
type = ADElementAverageSecondTimeDerivative
variable = v
[]
[]
[Outputs]
csv = true
[]
(test/tests/fvkernels/constraints/bounded_value.i)
[Mesh]
type = GeneratedMesh
dim = 1
nx = 20
[]
[Variables]
[v]
type = MooseVariableFVReal
# breaks the constraint
initial_condition = -1
[]
[lambda]
type = MooseVariableScalar
[]
[]
[FVKernels]
[time]
type = FVTimeKernel
variable = v
[]
[diff]
type = FVDiffusion
variable = v
coeff = coeff
[]
[average]
type = FVBoundedValueConstraint
variable = v
phi0 = 0
lambda = lambda
bound_type = 'HIGHER_THAN'
[]
[]
[FVBCs]
[left]
type = FVDirichletBC
variable = v
boundary = left
value = 7
[]
[right]
type = FVDirichletBC
variable = v
boundary = right
value = 0
[]
[]
[Materials]
[diff]
type = ADGenericFunctorMaterial
prop_names = 'coeff'
prop_values = '1'
[]
[]
[Executioner]
type = Transient
solve_type = 'Newton'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
num_steps = 2
dt = 0.001
nl_abs_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(test/tests/outputs/debug/show_execution_fv_flux_objects.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 = '2'
[]
[middle_boundary]
type = SideSetsBetweenSubdomainsGenerator
input = 'left'
primary_block = '0'
paired_block = '2'
new_boundary = 'middle'
[]
[]
[Variables]
[v]
family = MONOMIAL
order = CONSTANT
fv = true
block = 2
[]
[u]
type = MooseVariableFVReal
[]
[]
[ICs]
[v_ic]
type = FunctionIC
variable = v
function = 'if (x > 2 & x < 3, 0.5, 0)'
[]
[]
[FVKernels]
# Twice the kernel makes it not the Burgers equation, but shows the ordering
[2_burger]
type = FVBurgers1D
variable = v
[]
[1_burgers]
type = FVBurgers1D
variable = v
[]
[time]
type = FVTimeKernel
variable = v
[]
[time_u]
type = FVTimeKernel
variable = u
[]
[]
[FVBCs]
[fv_burgers_right]
type = FVBurgersOutflowBC
variable = v
boundary = 'middle'
[]
[fv_burgers_left]
type = FVBurgersOutflowBC
variable = v
boundary = 'left'
[]
[]
[FVInterfaceKernels]
[diff_ik]
type = FVOnlyAddDiffusionToOneSideOfInterface
variable1 = 'v'
variable2 = 'u'
boundary = 'middle'
coeff2 = '1'
subdomain1 = '2'
subdomain2 = '0'
[]
[]
[Executioner]
type = Transient
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
petsc_options = '-snes_converged_reason'
nl_abs_tol = 1e-7
nl_rel_tol = 1e-8
num_steps = 1
dt = 0.05
nl_forced_its = 1
line_search = none
[]
[Debug]
show_execution_order = 'LINEAR'
[]
(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
[]
(test/tests/indicators/gradient_jump_indicator/gradient_jump_indicator_fv_test.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 2
xmax = 2
nx = 2
ny = 1
subdomain_ids = '0 1'
[]
[interface_mesh]
type = SideSetsBetweenSubdomainsGenerator
input = gmg
primary_block = 0
paired_block = 1
new_boundary = interface
[]
# This creates enough elements to have defined gradients
[refine]
type = RefineBlockGenerator
input = interface_mesh
block = '0 1'
refinement = '3 3'
[]
[]
[Adaptivity]
marker = error_frac
max_h_level = 5
[Indicators]
[u0_jump]
type = GradientJumpIndicator
variable = u0
scale_by_flux_faces = false
[]
[]
[Markers]
[error_frac]
type = ErrorFractionMarker
coarsen = 0.15
indicator = u0_jump
refine = 0.7
[]
[]
[]
[Variables]
[u0]
family = MONOMIAL
order = CONSTANT
fv = true
block = 0
initial_condition = 0
[]
[u1]
family = MONOMIAL
order = CONSTANT
fv = true
block = 1
initial_condition = 0
[]
[]
[FVKernels]
[time0]
type = FVTimeKernel
variable = u0
[]
[diff0]
type = FVDiffusion
variable = u0
coeff = 1
block = 0
[]
[time1]
type = FVTimeKernel
variable = u1
[]
[diff1]
type = FVDiffusion
variable = u1
coeff = 1
block = 1
[]
[]
[FVInterfaceKernels]
[diffusion]
type = FVDiffusionInterface
variable1 = u0
variable2 = u1
boundary = interface
subdomain1 = 0
subdomain2 = 1
coeff1 = 1
coeff2 = 1
[]
[]
[FVBCs]
[left] # arbitrary user-chosen name
type = FVDirichletBC
variable = u0
boundary = 'left' # This must match a named boundary in the mesh file
value = 1
[]
[right] # arbitrary user-chosen name
type = FVNeumannBC
variable = u1
boundary = 'right' # This must match a named boundary in the mesh file
value = 0
[]
[]
[Executioner]
type = Transient
solve_type = 'Newton'
end_time = 0.5
dt = 0.1
[]
[VectorPostprocessors]
[samples]
type = LineValueSampler
variable = u0
# Avoiding element faces
start_point = '0.0001 1e-6 0'
end_point = '0.999999 1e-6 0'
num_points = 10
sort_by = 'x'
[]
[]
[Outputs]
execute_on = 'timestep_end'
csv = true
[]
(test/tests/fvkernels/fv_simple_diffusion/grad-adaptive.i)
[Mesh]
type = GeneratedMesh
dim = 1
nx = 100
[]
[Variables]
[v]
type = MooseVariableFVReal
initial_condition = 0
[]
[]
[AuxVariables]
[dummy]
type = MooseVariableFVReal
[]
[]
[FVKernels]
[time]
type = FVTimeKernel
variable = v
[]
[diff]
type = FVDiffusion
variable = v
coeff = coeff
[]
[]
[FVBCs]
[left]
type = FVDirichletBC
variable = v
boundary = left
value = 0
[]
[right]
type = FVDirichletBC
variable = v
boundary = right
value = 1
[]
[]
[Materials]
[diff]
type = ADGenericFunctorMaterial
prop_names = 'coeff'
prop_values = '1'
[]
[]
[Postprocessors]
[average]
type = ElementAverageValue
variable = v
[]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[TimeStepper]
type = IterationAdaptiveDT
dt = 1e-6
optimal_iterations = 6
[]
end_time = 1000
nl_abs_tol = 1e-8
[]
[Outputs]
exodus = false
[csv]
type = CSV
execute_on = 'final'
[]
[]
[Adaptivity]
steps = 1
marker = error
[Indicators]
[jump]
type = GradientJumpIndicator
variable = v
[]
[]
[Markers]
[error]
type = ErrorFractionMarker
coarsen = 0.1
refine = 0.7
indicator = jump
[]
[]
max_h_level = 1
[]
(modules/navier_stokes/test/tests/finite_volume/cns/benchmark_shock_tube_1D/hllc_sod_shocktube.i)
rho_left = 1
E_left = 2.501505578
u_left = 1e-15
rho_right = 0.125
E_right = 1.999770935
u_right = 1e-15
middle = 50
[GlobalParams]
fp = fp
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 1
xmin = 0
xmax = ${fparse 2 * middle}
nx = 1000
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[Variables]
[rho]
order = CONSTANT
family = MONOMIAL
fv = true
[]
[rho_u]
order = CONSTANT
family = MONOMIAL
fv = true
[]
[rho_E]
order = CONSTANT
family = MONOMIAL
fv = true
[]
[]
[AuxVariables]
[rho_a]
order = CONSTANT
family = MONOMIAL
[]
[]
[FVKernels]
[mass_time]
type = FVTimeKernel
variable = rho
[]
[mass_advection]
type = CNSFVMassHLLC
variable = rho
[]
[momentum_time]
type = FVTimeKernel
variable = rho_u
[]
[momentum_advection]
type = CNSFVMomentumHLLC
variable = rho_u
momentum_component = x
[]
[fluid_energy_time]
type = FVTimeKernel
variable = rho_E
[../]
[fluid_energy_advection]
type = CNSFVFluidEnergyHLLC
variable = rho_E
[]
[]
[FVBCs]
[mass_implicit]
type = CNSFVHLLCMassImplicitBC
variable = rho
fp = fp
boundary = 'left right'
[]
[mom_implicit]
type = CNSFVHLLCMomentumImplicitBC
variable = rho_u
momentum_component = x
fp = fp
boundary = 'left right'
[]
[fluid_energy_implicit]
type = CNSFVHLLCFluidEnergyImplicitBC
variable = rho_E
fp = fp
boundary = 'left right'
[]
[]
[ICs]
[rho_ic]
type = FunctionIC
variable = rho
function = 'if (x < ${middle}, ${rho_left}, ${rho_right})'
[]
[rho_u_ic]
type = FunctionIC
variable = rho_u
function = 'if (x < ${middle}, ${fparse rho_left * u_left}, ${fparse rho_right * u_right})'
[]
[rho_E_ic]
type = FunctionIC
variable = rho_E
function = 'if (x < ${middle}, ${fparse E_left * rho_left}, ${fparse E_right * rho_right})'
[]
[]
[Materials]
[var_mat]
type = ConservedVarValuesMaterial
rho = rho
rhou = rho_u
rho_et = rho_E
fp = fp
[]
[]
[Preconditioning]
active = ''
[./smp]
type = SMP
full = true
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
[../]
[]
[Executioner]
type = Transient
[TimeIntegrator]
type = ExplicitSSPRungeKutta
order = 2
[]
l_tol = 1e-8
start_time = 0.0
dt = 1e-2
end_time = 20
abort_on_solve_fail = true
[]
[Outputs]
exodus = true
perf_graph = true
[]
(test/tests/postprocessors/side_advection_flux_integral/side_advection_flux_integral_fv.i)
[Mesh]
[cmg]
type = CartesianMeshGenerator
dim = 2
dx = '0.75 0.75 0.75'
dy = '0.75 0.75 0.75'
ix = '2 2 2'
iy = '2 2 2'
subdomain_id = '1 1 1
1 2 1
1 1 1'
[]
[add_inner_boundaries_top]
type = SideSetsAroundSubdomainGenerator
input = cmg
new_boundary = 'block_2_top'
block = 2
normal = '0 1 0'
[]
[add_inner_boundaries_bot]
type = SideSetsAroundSubdomainGenerator
input = add_inner_boundaries_top
new_boundary = 'block_2_bot'
block = 2
normal = '0 -1 0'
[]
[add_inner_boundaries_right]
type = SideSetsAroundSubdomainGenerator
input = add_inner_boundaries_bot
new_boundary = 'block_2_right'
block = 2
normal = '1 0 0'
[]
[add_inner_boundaries_left]
type = SideSetsAroundSubdomainGenerator
input = add_inner_boundaries_right
new_boundary = 'block_2_left'
block = 2
normal = '-1 0 0'
[]
[]
[Variables]
[u]
type = MooseVariableFVReal
two_term_boundary_expansion = false
[]
[]
[ICs]
[u_blob]
type = FunctionIC
variable = u
function = 'if(x<0.75,if(y<0.75,1,0),0)'
[]
[]
[FVKernels]
[advection]
type = FVAdvection
variable = u
velocity = '2 1.5 0'
[]
[time]
type = FVTimeKernel
variable = u
[]
[]
[FVBCs]
[fv_outflow]
type = FVConstantScalarOutflowBC
velocity = '2 1.5 0'
variable = u
boundary = 'right top'
[]
[]
[Executioner]
type = Transient
solve_type = LINEAR
dt = 0.01
end_time = 0.02
l_tol = 1E-14
[]
[Postprocessors]
[flux_right]
type = SideAdvectiveFluxIntegral
boundary = 'block_2_right'
vel_x = 2
vel_y = 1.5
component = x
advected_quantity = u
[]
[flux_left_exact]
type = SideAdvectiveFluxIntegral
boundary = 'block_2_left'
vel_x = 2
vel_y = 1.5
component = x
advected_quantity = u
[]
[]
[Outputs]
csv = true
[]
(test/tests/postprocessors/pseudotimestep/fv_burgers_pseudo.i)
[Mesh]
[gen_mesh]
type = GeneratedMeshGenerator
dim = 1
xmin = -1
xmax = 1
nx = 500
[]
[]
[Variables]
[v]
family = MONOMIAL
order = CONSTANT
fv = true
[]
[]
[ICs]
[v_ic]
type = FunctionIC
variable = v
function = '-1/(1+exp(-(x-z)/2/0.0005))'
[]
[]
[FVKernels]
[burgers]
type = FVBurgers1D
variable = v
[]
[difussion]
type = FVDiffusion
coeff= 0.0005
variable = v
[]
[time]
type = FVTimeKernel
variable = v
[]
[]
[FVBCs]
[fv_burgers_outflow]
type = FVBurgersOutflowBC
variable = v
boundary = 'left right'
[]
[]
[Postprocessors]
[pseudotimestep]
type = PseudoTimestep
method = 'SER'
initial_dt = 1
alpha = 1.5
[]
[]
[Executioner]
type = Transient
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
petsc_options = '-snes_converged_reason'
num_steps = 5
[TimeStepper]
type = PostprocessorDT
postprocessor = pseudotimestep
[]
[]
[Outputs]
print_linear_residuals = false
csv = true
[]
(test/tests/auxkernels/advection_flux/advection_flux_fv.i)
[Mesh]
type = GeneratedMesh # Can generate simple lines, rectangles and rectangular prisms
dim = 2 # Dimension of the mesh
nx = 10 # Number of elements in the x direction
ny = 10 # Number of elements in the y direction
xmax = 1.0
ymax = 1.0
[]
[Variables]
[u]
type = MooseVariableFVReal
two_term_boundary_expansion = false
[]
[]
[AuxVariables]
[flux_x]
type = MooseVariableFVReal
order = CONSTANT
family = MONOMIAL
[]
[]
[ICs]
[u_ic]
type = FunctionIC
variable = u
function = 'r2 := (x - 0.5)*(x - 0.5) + (y - 0.3)*(y - 0.3); exp(-r2 * 20)'
[]
[]
[FVKernels]
[advection]
type = FVAdvection
variable = u
velocity = '1 0.5 0'
[]
[time]
type = FVTimeKernel
variable = u
[]
[]
[FVBCs]
[fv_outflow]
type = FVConstantScalarOutflowBC
velocity = '1 0.5 0'
variable = u
boundary = 'right top'
[]
[]
[AuxKernels]
[flux_x]
type = AdvectiveFluxAux
variable = flux_x
vel_x = 1
vel_y = 0.5
advected_variable = u
component = normal
boundary = 'left right'
check_boundary_restricted = false
[]
[]
[Postprocessors]
[flux_right]
type = SideIntegralVariablePostprocessor
variable = flux_x
boundary = 'right'
[]
[flux_right_exact]
type = SideAdvectiveFluxIntegral
vel_x = 1
vel_y = 0.5
component = normal
advected_quantity = u
boundary = 'right'
[]
[flux_left]
type = SideIntegralVariablePostprocessor
variable = flux_x
boundary = 'left'
[]
[flux_left_exact]
type = SideAdvectiveFluxIntegral
vel_x = 1
vel_y = 0.5
component = normal
advected_quantity = u
boundary = 'left'
[]
[]
[Executioner]
type = Transient
petsc_options = '-snes_converged_reason'
num_steps = 10
dt = 0.25
nl_rel_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(test/tests/fvkernels/fv_constant_scalar_advection/2D_constant_scalar_advection.i)
[Mesh]
[gen_mesh]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 2
ymin = 0
ymax = 4
nx = 10
ny = 20
[]
[]
[Variables]
[v]
type = MooseVariableFVReal
two_term_boundary_expansion = false
[]
[]
[ICs]
[v_ic]
type = FunctionIC
variable = v
function = 'r2 := (x - 0.5)*(x - 0.5) + (y - 0.3)*(y - 0.3); exp(-r2 * 20)'
[]
[]
[FVKernels]
[advection]
type = FVAdvection
variable = v
velocity = '1 0.5 0'
[]
[time]
type = FVTimeKernel
variable = v
[]
[]
[FVBCs]
[fv_outflow]
type = FVConstantScalarOutflowBC
velocity = '1 0.5 0'
variable = v
boundary = 'right top'
[]
[]
[Executioner]
type = Transient
petsc_options = '-snes_converged_reason'
num_steps = 10
dt = 0.25
nl_rel_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/cns/implicit_bcs/hllc_sod_shocktube.i)
rho_left = 1
E_left = 2.501505578
u_left = 1e-15
rho_right = 0.125
E_right = 1.999770935
u_right = 1e-15
middle = 0.5
[GlobalParams]
fp = fp
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${fparse 2 * middle}
nx = 5
ymin = 0
ymax = 1
ny = 2
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
allow_imperfect_jacobians = true
[]
[]
[Variables]
[rho]
order = CONSTANT
family = MONOMIAL
fv = true
[]
[rho_u]
order = CONSTANT
family = MONOMIAL
fv = true
[]
[rho_v]
order = CONSTANT
family = MONOMIAL
fv = true
initial_condition = 1e-10
[]
[rho_E]
order = CONSTANT
family = MONOMIAL
fv = true
[]
[]
[FVKernels]
[mass_time]
type = FVTimeKernel
variable = rho
[]
[mass_advection]
type = CNSFVMassHLLC
variable = rho
[]
[momentum_x_time]
type = FVTimeKernel
variable = rho_u
[]
[momentum_x_advection]
type = CNSFVMomentumHLLC
variable = rho_u
momentum_component = x
[]
[momentum_y_time]
type = FVTimeKernel
variable = rho_v
[]
[momentum_y_advection]
type = CNSFVMomentumHLLC
variable = rho_v
momentum_component = y
[]
[fluid_energy_time]
type = FVTimeKernel
variable = rho_E
[]
[fluid_energy_advection]
type = CNSFVFluidEnergyHLLC
variable = rho_E
[]
[]
[FVBCs]
[mass_implicit]
type = CNSFVHLLCMassImplicitBC
variable = rho
fp = fp
boundary = 'left right'
[]
[mom_x_implicit]
type = CNSFVHLLCMomentumImplicitBC
variable = rho_u
momentum_component = x
fp = fp
boundary = 'left right'
[]
[wall]
type = CNSFVMomImplicitPressureBC
variable = rho_v
momentum_component = y
boundary = 'top bottom'
[]
[fluid_energy_implicit]
type = CNSFVHLLCFluidEnergyImplicitBC
variable = rho_E
fp = fp
boundary = 'left right'
[]
[]
[ICs]
[rho_ic]
type = FunctionIC
variable = rho
function = 'if (x < ${middle}, ${rho_left}, ${rho_right})'
[]
[rho_u_ic]
type = FunctionIC
variable = rho_u
function = 'if (x < ${middle}, ${fparse rho_left * u_left}, ${fparse rho_right * u_right})'
[]
[rho_E_ic]
type = FunctionIC
variable = rho_E
function = 'if (x < ${middle}, ${fparse E_left * rho_left}, ${fparse E_right * rho_right})'
[]
[]
[Materials]
[var_mat]
type = ConservedVarValuesMaterial
rho = rho
rhou = rho_u
rhov = rho_v
rho_et = rho_E
fp = fp
[]
[]
[Executioner]
type = Transient
[TimeIntegrator]
type = ExplicitSSPRungeKutta
order = 2
[]
l_tol = 1e-8
# run to t = 0.15
start_time = 0.0
dt = 1e-1
end_time = 10
abort_on_solve_fail = true
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/cns/stagnation_inlet/supersonic_nozzle_hllc.i)
stagnation_pressure = 1
stagnation_temperature = 1
[GlobalParams]
fp = fp
[]
[Debug]
show_material_props = true
[]
[Mesh]
[file]
type = FileMeshGenerator
file = supersonic_nozzle.e
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[Variables]
[rho]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 0.0034
[]
[rho_u]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 1e-4
outputs = none
[]
[rho_v]
family = MONOMIAL
order = CONSTANT
fv = true
outputs = none
[]
[rho_E]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 2.5
[]
[]
[FVKernels]
# Mass conservation
[mass_time]
type = FVTimeKernel
variable = rho
[]
[mass_advection]
type = CNSFVMassHLLC
variable = rho
[]
# Momentum x conservation
[momentum_x_time]
type = FVTimeKernel
variable = rho_u
[]
[momentum_x_advection]
type = CNSFVMomentumHLLC
variable = rho_u
momentum_component = x
[]
# Momentum y conservation
[momentum_y_time]
type = FVTimeKernel
variable = rho_v
[]
[momentum_y_advection]
type = CNSFVMomentumHLLC
variable = rho_v
momentum_component = y
[]
# Fluid energy conservation
[fluid_energy_time]
type = FVTimeKernel
variable = rho_E
[]
[fluid_energy_advection]
type = CNSFVFluidEnergyHLLC
variable = rho_E
[]
[]
[FVBCs]
## inflow stagnation boundaries
[mass_stagnation_inflow]
type = CNSFVHLLCMassStagnationInletBC
variable = rho
stagnation_pressure = ${stagnation_pressure}
stagnation_temperature = ${stagnation_temperature}
boundary = left
[]
[momentum_x_stagnation_inflow]
type = CNSFVHLLCMomentumStagnationInletBC
variable = rho_u
momentum_component = x
stagnation_pressure = ${stagnation_pressure}
stagnation_temperature = ${stagnation_temperature}
boundary = left
[]
[momentum_y_stagnation_inflow]
type = CNSFVHLLCMomentumStagnationInletBC
variable = rho_v
momentum_component = y
stagnation_pressure = ${stagnation_pressure}
stagnation_temperature = ${stagnation_temperature}
boundary = left
[../]
[fluid_energy_stagnation_inflow]
type = CNSFVHLLCFluidEnergyStagnationInletBC
variable = rho_E
stagnation_pressure = ${stagnation_pressure}
stagnation_temperature = ${stagnation_temperature}
boundary = left
[]
## outflow implicit conditions
[mass_outflow]
type = CNSFVHLLCMassImplicitBC
variable = rho
boundary = right
[]
[momentum_x_outflow]
type = CNSFVHLLCMomentumImplicitBC
variable = rho_u
momentum_component = x
boundary = right
[]
[momentum_y_outflow]
type = CNSFVHLLCMomentumImplicitBC
variable = rho_v
momentum_component = y
boundary = right
[]
[fluid_energy_outflow]
type = CNSFVHLLCFluidEnergyImplicitBC
variable = rho_E
boundary = right
[]
# wall conditions
[momentum_x_pressure_wall]
type = CNSFVMomImplicitPressureBC
variable = rho_u
momentum_component = x
boundary = wall
[]
[momentum_y_pressure_wall]
type = CNSFVMomImplicitPressureBC
variable = rho_v
momentum_component = y
boundary = wall
[]
[]
[AuxVariables]
[Ma]
family = MONOMIAL
order = CONSTANT
[]
[Ma_layered]
family = MONOMIAL
order = CONSTANT
[]
[]
[UserObjects]
[layered_Ma_UO]
type = LayeredAverage
variable = Ma
num_layers = 100
direction = x
[]
[]
[AuxKernels]
[Ma_aux]
type = NSMachAux
variable = Ma
fluid_properties = fp
use_material_properties = true
[]
[Ma_layered_aux]
type = SpatialUserObjectAux
variable = Ma_layered
user_object = layered_Ma_UO
[]
[]
[Materials]
[var_mat]
type = ConservedVarValuesMaterial
rho = rho
rhou = rho_u
rhov = rho_v
rho_et = rho_E
[]
[fluid_props]
type = GeneralFluidProps
porosity = 1
characteristic_length = 1
[]
[sound_speed]
type = SoundspeedMat
fp = fp
[]
[]
[Postprocessors]
[cfl_dt]
type = ADCFLTimeStepSize
c_names = 'sound_speed'
vel_names = 'speed'
CFL = 0.5
[]
[outflow_Ma]
type = SideAverageValue
variable = Ma
boundary = right
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
[]
[]
[Executioner]
type = Transient
end_time = 0.1
[TimeIntegrator]
type = ExplicitSSPRungeKutta
order = 2
[]
l_tol = 1e-8
[TimeStepper]
type = PostprocessorDT
postprocessor = cfl_dt
[]
[]
[VectorPostprocessors]
[Ma_layered]
type = LineValueSampler
variable = Ma_layered
start_point = '0 0 0'
end_point = '10 0 0'
num_points = 100
sort_by = x
[]
[]
[Outputs]
exodus = true
[]
(modules/porous_flow/test/tests/heat_conduction/two_phase_fv.i)
# 2 phase heat conduction, with saturation fixed at 0.5
# Apply a boundary condition of T=300 to a bar that
# is initially at T=200, and observe the expected
# error-function response
[Mesh]
[mesh]
type = GeneratedMeshGenerator
dim = 1
nx = 10
xmin = 0
xmax = 100
[]
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[phase0_porepressure]
type = MooseVariableFVReal
initial_condition = 0
[]
[phase1_saturation]
type = MooseVariableFVReal
initial_condition = 0.5
[]
[temp]
type = MooseVariableFVReal
initial_condition = 200
[]
[]
[FVKernels]
[dummy_p0]
type = FVTimeKernel
variable = phase0_porepressure
[]
[dummy_s1]
type = FVTimeKernel
variable = phase1_saturation
[]
[energy_dot]
type = FVPorousFlowEnergyTimeDerivative
variable = temp
[]
[heat_conduction]
type = FVPorousFlowHeatConduction
variable = temp
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'temp phase0_porepressure phase1_saturation'
number_fluid_phases = 2
number_fluid_components = 1
[]
[pc]
type = PorousFlowCapillaryPressureConst
pc = 0
[]
[]
[FluidProperties]
[simple_fluid0]
type = SimpleFluidProperties
bulk_modulus = 1.5
density0 = 0.4
thermal_expansion = 0
cv = 1
[]
[simple_fluid1]
type = SimpleFluidProperties
bulk_modulus = 0.5
density0 = 0.3
thermal_expansion = 0
cv = 2
[]
[]
[Materials]
[temperature]
type = ADPorousFlowTemperature
temperature = temp
[]
[thermal_conductivity]
type = ADPorousFlowThermalConductivityIdeal
dry_thermal_conductivity = '0.3 0 0 0 0 0 0 0 0'
wet_thermal_conductivity = '1.7 0 0 0 0 0 0 0 0'
exponent = 1.0
aqueous_phase_number = 1
[]
[ppss]
type = ADPorousFlow2PhasePS
phase0_porepressure = phase0_porepressure
phase1_saturation = phase1_saturation
capillary_pressure = pc
[]
[porosity]
type = ADPorousFlowPorosityConst
porosity = 0.8
[]
[rock_heat]
type = ADPorousFlowMatrixInternalEnergy
specific_heat_capacity = 1.0
density = 0.25
[]
[simple_fluid0]
type = ADPorousFlowSingleComponentFluid
fp = simple_fluid0
phase = 0
[]
[simple_fluid1]
type = ADPorousFlowSingleComponentFluid
fp = simple_fluid1
phase = 1
[]
[]
[FVBCs]
[left]
type = FVDirichletBC
boundary = left
value = 300
variable = temp
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 1E1
end_time = 1E2
[]
[Postprocessors]
[t005]
type = PointValue
variable = temp
point = '5 0 0'
execute_on = 'initial timestep_end'
[]
[t015]
type = PointValue
variable = temp
point = '15 0 0'
execute_on = 'initial timestep_end'
[]
[t025]
type = PointValue
variable = temp
point = '25 0 0'
execute_on = 'initial timestep_end'
[]
[t035]
type = PointValue
variable = temp
point = '35 0 0'
execute_on = 'initial timestep_end'
[]
[t045]
type = PointValue
variable = temp
point = '45 0 0'
execute_on = 'initial timestep_end'
[]
[t055]
type = PointValue
variable = temp
point = '55 0 0'
execute_on = 'initial timestep_end'
[]
[t065]
type = PointValue
variable = temp
point = '65 0 0'
execute_on = 'initial timestep_end'
[]
[t075]
type = PointValue
variable = temp
point = '75 0 0'
execute_on = 'initial timestep_end'
[]
[t085]
type = PointValue
variable = temp
point = '85 0 0'
execute_on = 'initial timestep_end'
[]
[t095]
type = PointValue
variable = temp
point = '95 0 0'
execute_on = 'initial timestep_end'
[]
[]
[Outputs]
file_base = two_phase_fv
csv = true
[]
(modules/navier_stokes/test/tests/finite_volume/cns/shock_tube_2D_cavity/hllc_sod_shocktube_2D.i)
rho_left = 1
E_left = 2.501505578
u_left = 1e-15
rho_right = 0.125
E_right = 1.999770935
u_right = 1e-15
x_sep = 35
[GlobalParams]
fp = fp
[]
[Mesh]
[./cartesian]
type = CartesianMeshGenerator
dim = 2
dx = '40 20'
ix = '200 100'
dy = '1 20 2 20 1'
iy = '4 100 10 100 4'
subdomain_id = '0 0
0 1
1 1
0 1
0 0'
[../]
[./wall]
type = SideSetsBetweenSubdomainsGenerator
input = cartesian
primary_block = 1
paired_block = 0
new_boundary = 'wall'
[../]
[./delete]
type = BlockDeletionGenerator
input = wall
block = 0
[../]
[]
[FluidProperties]
[./fp]
type = IdealGasFluidProperties
allow_imperfect_jacobians = true
[../]
[]
[Variables]
[./rho]
order = CONSTANT
family = MONOMIAL
fv = true
[../]
[./rho_u]
order = CONSTANT
family = MONOMIAL
fv = true
[../]
[./rho_v]
order = CONSTANT
family = MONOMIAL
fv = true
[../]
[./rho_E]
order = CONSTANT
family = MONOMIAL
fv = true
[../]
[]
[AuxVariables]
[./Ma]
order = CONSTANT
family = MONOMIAL
[../]
[./p]
order = CONSTANT
family = MONOMIAL
[../]
[./v_norm]
order = CONSTANT
family = MONOMIAL
[../]
[./temperature]
order = CONSTANT
family = MONOMIAL
[../]
[]
[AuxKernels]
[./Ma_aux]
type = NSMachAux
variable = Ma
fluid_properties = fp
use_material_properties = true
[../]
[./p_aux]
type = ADMaterialRealAux
variable = p
property = pressure
[../]
[./v_norm_aux]
type = ADMaterialRealAux
variable = v_norm
property = speed
[../]
[./temperature_aux]
type = ADMaterialRealAux
variable = temperature
property = T_fluid
[../]
[]
[FVKernels]
[./mass_time]
type = FVTimeKernel
variable = rho
[../]
[./mass_advection]
type = CNSFVMassHLLC
variable = rho
[../]
[./momentum_x_time]
type = FVTimeKernel
variable = rho_u
[../]
[./momentum_x_advection]
type = CNSFVMomentumHLLC
variable = rho_u
momentum_component = x
[../]
[./momentum_y_time]
type = FVTimeKernel
variable = rho_v
[../]
[./momentum_y_advection]
type = CNSFVMomentumHLLC
variable = rho_v
momentum_component = y
[../]
[./fluid_energy_time]
type = FVTimeKernel
variable = rho_E
[../]
[./fluid_energy_advection]
type = CNSFVFluidEnergyHLLC
variable = rho_E
[../]
[]
[FVBCs]
[./mom_x_pressure]
type = CNSFVMomImplicitPressureBC
variable = rho_u
momentum_component = x
boundary = 'left right wall'
[../]
[./mom_y_pressure]
type = CNSFVMomImplicitPressureBC
variable = rho_v
momentum_component = y
boundary = 'wall'
[../]
[]
[ICs]
[./rho_ic]
type = FunctionIC
variable = rho
function = 'if (x < ${x_sep}, ${rho_left}, ${rho_right})'
[../]
[./rho_u_ic]
type = FunctionIC
variable = rho_u
function = 'if (x < ${x_sep}, ${fparse rho_left * u_left}, ${fparse rho_right * u_right})'
[../]
[./rho_E_ic]
type = FunctionIC
variable = rho_E
function = 'if (x < ${x_sep}, ${fparse E_left * rho_left}, ${fparse E_right * rho_right})'
[../]
[]
[Materials]
[./var_mat]
type = ConservedVarValuesMaterial
rho = rho
rhou = rho_u
rhov = rho_v
rho_et = rho_E
fp = fp
[../]
[./sound_speed]
type = SoundspeedMat
fp = fp
[../]
[]
[Preconditioning]
[./smp]
type = SMP
full = true
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
[../]
[]
[Postprocessors]
[./cfl_dt]
type = ADCFLTimeStepSize
c_names = 'sound_speed'
vel_names = 'speed'
[../]
[]
[Executioner]
type = Transient
end_time = 100
[TimeIntegrator]
type = ExplicitSSPRungeKutta
order = 2
[]
l_tol = 1e-8
[./TimeStepper]
type = PostprocessorDT
postprocessor = cfl_dt
[../]
[]
(test/tests/tag/tag-fv.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 5
ny = 5
[]
[Variables]
[v]
type = MooseVariableFVReal
initial_condition = 7
[]
[]
[AuxVariables]
[soln_dof]
type = MooseVariableFVReal
[]
[soln_old_dof]
type = MooseVariableFVReal
[]
[soln_older_dof]
type = MooseVariableFVReal
[]
[resid_nontime_dof]
type = MooseVariableFVReal
[]
[soln]
type = MooseVariableFVReal
[]
[soln_old]
type = MooseVariableFVReal
[]
[soln_older]
type = MooseVariableFVReal
[]
[resid_nontime]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[soln_dof]
type = TagVectorDofValueAux
variable = soln_dof
v = v
vector_tag = 'solution'
[]
[soln_old_dof]
type = TagVectorDofValueAux
variable = soln_old_dof
v = v
vector_tag = 'solution_state_1'
[]
[soln_older_dof]
type = TagVectorDofValueAux
variable = soln_older_dof
v = v
vector_tag = 'solution_state_2'
[]
[nontime_dof]
type = TagVectorDofValueAux
variable = resid_nontime_dof
v = v
vector_tag = 'nontime'
[]
[soln]
type = TagVectorAux
variable = soln
v = v
vector_tag = 'solution'
[]
[soln_old]
type = TagVectorAux
variable = soln_old
v = v
vector_tag = 'solution_state_1'
[]
[soln_older]
type = TagVectorAux
variable = soln_older
v = v
vector_tag = 'solution_state_2'
[]
[nontime]
type = TagVectorAux
variable = resid_nontime
v = v
vector_tag = 'nontime'
[]
[]
[FVKernels]
[time]
type = FVTimeKernel
variable = v
[]
[diff]
type = FVDiffusion
variable = v
coeff = coeff
[]
[]
[FVBCs]
[left]
type = FVDirichletBC
variable = v
boundary = left
value = 7
[]
[right]
type = FVDirichletBC
variable = v
boundary = right
value = 42
[]
[]
[Materials]
[diff]
type = ADGenericFunctorMaterial
prop_names = 'coeff'
prop_values = '.2'
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
num_steps = 5
dt = 0.1
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/cns/pressure_outlet/subsonic_nozzle_fixed_inflow_hllc.i)
inlet_vel = 120
rho_in = 0.8719748696
H_in = 4.0138771448e+05
gamma = 1.4
R = 8.3145
molar_mass = 29e-3
R_specific = ${fparse R / molar_mass}
cp = ${fparse gamma * R_specific / (gamma - 1)}
cv = ${fparse cp / gamma}
T_in = ${fparse H_in / gamma / cv}
mass_flux = ${fparse inlet_vel * rho_in}
outlet_pressure = 0.9e5
[GlobalParams]
fp = fp
[]
[Debug]
show_material_props = true
[]
[Mesh]
[file]
type = FileMeshGenerator
file = subsonic_nozzle.e
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[Variables]
[rho]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 0.8719748696
[]
[rho_u]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 1e-4
[]
[rho_v]
family = MONOMIAL
order = CONSTANT
fv = true
[]
[rho_E]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 2.5e5
[]
[]
[FVKernels]
# Mass conservation
[mass_time]
type = FVTimeKernel
variable = rho
[]
[mass_advection]
type = CNSFVMassHLLC
variable = rho
[]
# Momentum x conservation
[momentum_x_time]
type = FVTimeKernel
variable = rho_u
[]
[momentum_x_advection]
type = CNSFVMomentumHLLC
variable = rho_u
momentum_component = x
[]
# Momentum y conservation
[momentum_y_time]
type = FVTimeKernel
variable = rho_v
[]
[momentum_y_advection]
type = CNSFVMomentumHLLC
variable = rho_v
momentum_component = y
[]
# Fluid energy conservation
[fluid_energy_time]
type = FVTimeKernel
variable = rho_E
[]
[fluid_energy_advection]
type = CNSFVFluidEnergyHLLC
variable = rho_E
[]
[]
[FVBCs]
## inflow boundaries
[mass_inflow]
type = CNSFVHLLCSpecifiedMassFluxAndTemperatureMassBC
variable = rho
boundary = left
rhou = ${mass_flux}
rhov = 0
temperature = ${T_in}
[]
[momentum_x_inflow]
type = CNSFVHLLCSpecifiedMassFluxAndTemperatureMomentumBC
variable = rho_u
boundary = left
rhou = ${mass_flux}
rhov = 0
temperature = ${T_in}
momentum_component = x
[]
[momentum_y_inflow]
type = CNSFVHLLCSpecifiedMassFluxAndTemperatureMomentumBC
variable = rho_v
boundary = left
rhou = ${mass_flux}
rhov = 0
temperature = ${T_in}
momentum_component = y
[]
[fluid_energy_inflow]
type = CNSFVHLLCSpecifiedMassFluxAndTemperatureFluidEnergyBC
variable = rho_E
boundary = left
rhou = ${mass_flux}
rhov = 0
temperature = ${T_in}
[]
## outflow conditions
[mass_outflow]
type = CNSFVHLLCSpecifiedPressureMassBC
variable = rho
boundary = right
pressure = ${outlet_pressure}
[]
[momentum_x_outflow]
type = CNSFVHLLCSpecifiedPressureMomentumBC
variable = rho_u
boundary = right
momentum_component = x
pressure = ${outlet_pressure}
[]
[momentum_y_outflow]
type = CNSFVHLLCSpecifiedPressureMomentumBC
variable = rho_v
boundary = right
momentum_component = y
pressure = ${outlet_pressure}
[]
[fluid_energy_outflow]
type = CNSFVHLLCSpecifiedPressureFluidEnergyBC
variable = rho_E
boundary = right
pressure = ${outlet_pressure}
[]
# wall conditions
[momentum_x_pressure_wall]
type = CNSFVMomImplicitPressureBC
variable = rho_u
momentum_component = x
boundary = wall
[]
[momentum_y_pressure_wall]
type = CNSFVMomImplicitPressureBC
variable = rho_v
momentum_component = y
boundary = wall
[]
[]
[AuxVariables]
[Ma]
family = MONOMIAL
order = CONSTANT
[]
[p]
family = MONOMIAL
order = CONSTANT
[]
[Ma_layered]
family = MONOMIAL
order = CONSTANT
[]
[]
[UserObjects]
[layered_Ma_UO]
type = LayeredAverage
variable = Ma
num_layers = 10
direction = x
[]
[]
[AuxKernels]
[Ma_aux]
type = NSMachAux
variable = Ma
fluid_properties = fp
use_material_properties = true
[]
[p_aux]
type = ADMaterialRealAux
variable = p
property = pressure
[]
[Ma_layered_aux]
type = SpatialUserObjectAux
variable = Ma_layered
user_object = layered_Ma_UO
[]
[]
[Materials]
[var_mat]
type = ConservedVarValuesMaterial
rho = rho
rhou = rho_u
rhov = rho_v
rho_et = rho_E
[]
[sound_speed]
type = SoundspeedMat
[]
[]
[Postprocessors]
[outflow_Ma]
type = SideAverageValue
variable = Ma
boundary = right
[]
[outflow_pressure]
type = SideAverageValue
variable = p
boundary = right
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
[]
[]
[Executioner]
type = Transient
end_time = 10
solve_type = NEWTON
nl_abs_tol = 1e-7
[TimeIntegrator]
type = ImplicitEuler
[]
[TimeStepper]
type = IterationAdaptiveDT
dt = 5e-3
optimal_iterations = 6
growth_factor = 1.5
[]
[]
[VectorPostprocessors]
[Ma_layered]
type = LineValueSampler
variable = Ma_layered
start_point = '0 0 0'
end_point = '3 0 0'
num_points = 100
sort_by = x
warn_discontinuous_face_values = false
[]
[]
[Outputs]
[out]
type = Exodus
execute_on = 'final'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/cns/symmetry_test/2D_symmetry.i)
rho_inside = 1
E_inside = 2.501505578
rho_outside = 0.125
E_outside = 1.999770935
radius = 0.1
angle = 45
[GlobalParams]
fp = fp
[]
[Debug]
show_material_props = true
[]
[Mesh]
[file]
type = GeneratedMeshGenerator
dim = 2
xmin = -0.5
xmax = 0.5
nx = 10
ymin = -0.5
ymax = 0.5
ny = 10
[../]
[rotate]
type = TransformGenerator
vector_value = '${angle} 0 0'
transform = ROTATE
input = file
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
allow_imperfect_jacobians = true
[]
[]
[Variables]
[rho]
family = MONOMIAL
order = CONSTANT
fv = true
[../]
[rho_u]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 1e-15
outputs = none
[]
[rho_v]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 1e-15
outputs = none
[]
[rho_E]
family = MONOMIAL
order = CONSTANT
fv = true
[]
[]
[ICs]
[rho_ic]
type = FunctionIC
variable = rho
function = 'if (abs(x) < ${radius} & abs(y) < ${radius}, ${rho_inside}, ${rho_outside})'
[]
[rho_E_ic]
type = FunctionIC
variable = rho_E
function = 'if (abs(x) < ${radius} & abs(y) < ${radius}, ${fparse E_inside * rho_inside}, ${fparse E_outside * rho_outside})'
[]
[]
[FVKernels]
# Mass conservation
[mass_time]
type = FVTimeKernel
variable = rho
[]
[mass_advection]
type = CNSFVMassHLLC
variable = rho
fp = fp
[]
# Momentum x conservation
[momentum_x_time]
type = FVTimeKernel
variable = rho_u
[]
[momentum_x_advection]
type = CNSFVMomentumHLLC
variable = rho_u
momentum_component = x
fp = fp
[]
# Momentum y conservation
[momentum_y_time]
type = FVTimeKernel
variable = rho_v
[]
[./momentum_y_advection]
type = CNSFVMomentumHLLC
variable = rho_v
momentum_component = y
[]
# Fluid energy conservation
[./fluid_energy_time]
type = FVTimeKernel
variable = rho_E
[]
[./fluid_energy_advection]
type = CNSFVFluidEnergyHLLC
variable = rho_E
fp = fp
[]
[]
[FVBCs]
## outflow implicit conditions
[mass_outflow]
type = CNSFVHLLCMassImplicitBC
variable = rho
fp = fp
boundary = 'left right top bottom'
[]
[./momentum_x_outflow]
type = CNSFVHLLCMomentumImplicitBC
variable = rho_u
momentum_component = x
fp = fp
boundary = 'left right top bottom'
[]
[momentum_y_outflow]
type = CNSFVHLLCMomentumImplicitBC
variable = rho_v
momentum_component = y
fp = fp
boundary = 'left right top bottom'
[]
[fluid_energy_outflow]
type = CNSFVHLLCFluidEnergyImplicitBC
variable = rho_E
fp = fp
boundary = 'left right top bottom'
[]
[]
[AuxVariables]
[Ma]
family = MONOMIAL
order = CONSTANT
[]
[p]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[Ma_aux]
type = NSMachAux
variable = Ma
fluid_properties = fp
use_material_properties = true
[]
[p_aux]
type = ADMaterialRealAux
variable = p
property = pressure
[]
[]
[Materials]
[var_mat]
type = ConservedVarValuesMaterial
rho = rho
rhou = rho_u
rhov = rho_v
rho_et = rho_E
[]
[sound_speed]
type = SoundspeedMat
fp = fp
[]
[]
[Postprocessors]
[cfl_dt]
type = ADCFLTimeStepSize
c_names = 'sound_speed'
vel_names = 'speed'
CFL = 0.5
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
[]
[]
[Executioner]
type = Transient
end_time = 0.2
[TimeIntegrator]
type = ExplicitSSPRungeKutta
order = 2
[]
l_tol = 1e-8
[TimeStepper]
type = PostprocessorDT
postprocessor = cfl_dt
[]
[]
(modules/porous_flow/include/fvkernels/FVPorousFlowMassTimeDerivative.h)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#pragma once
#include "FVTimeKernel.h"
class PorousFlowDictator;
/**
* Time derivative of fluid mass
*/
class FVPorousFlowMassTimeDerivative : public FVTimeKernel
{
public:
static InputParameters validParams();
FVPorousFlowMassTimeDerivative(const InputParameters & parameters);
protected:
ADReal computeQpResidual() override;
/// UserObject that holds information (number of phases, components, etc)
const PorousFlowDictator & _dictator;
/// Number of fluid phases present
const unsigned int _num_phases;
/// Index of the fluid component this kernel applies to
const unsigned int _fluid_component;
/// Porosity
const ADMaterialProperty<Real> & _porosity;
const MaterialProperty<Real> & _porosity_old;
/// Fluid density
const ADMaterialProperty<std::vector<Real>> & _density;
const MaterialProperty<std::vector<Real>> & _density_old;
/// Fluid phase saturation
const ADMaterialProperty<std::vector<Real>> & _saturation;
const MaterialProperty<std::vector<Real>> & _saturation_old;
/// Mass fraction of fluid components in fluid phases
const ADMaterialProperty<std::vector<std::vector<Real>>> & _mass_fractions;
const MaterialProperty<std::vector<std::vector<Real>>> & _mass_fractions_old;
};
(modules/navier_stokes/include/fvkernels/FVPorosityTimeDerivative.h)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#pragma once
#include "FVTimeKernel.h"
/**
* Applies a residual equal to
* \f$\epsilon \frac{\partial u}{\partial t}\f$
*/
class FVPorosityTimeDerivative : public FVTimeKernel
{
public:
static InputParameters validParams();
FVPorosityTimeDerivative(const InputParameters & parameters);
protected:
ADReal computeQpResidual() override;
/// The porosity
const MaterialProperty<Real> & _eps;
};
(modules/navier_stokes/include/fvkernels/PCNSFVDensityTimeDerivative.h)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#pragma once
#include "FVTimeKernel.h"
class PCNSFVDensityTimeDerivative : public FVTimeKernel
{
public:
static InputParameters validParams();
PCNSFVDensityTimeDerivative(const InputParameters & parameters);
protected:
ADReal computeQpResidual() override;
/// The time derivative of the primary variable
const ADVariableValue & _u_dot;
/// The porosity
const MaterialProperty<Real> & _eps;
const ADVariableValue & _rho_dot;
const ADVariableValue & _rho;
};
(modules/porous_flow/include/fvkernels/FVPorousFlowEnergyTimeDerivative.h)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#pragma once
#include "FVTimeKernel.h"
class PorousFlowDictator;
/**
* Time derivative of energy
*/
class FVPorousFlowEnergyTimeDerivative : public FVTimeKernel
{
public:
static InputParameters validParams();
FVPorousFlowEnergyTimeDerivative(const InputParameters & parameters);
protected:
ADReal computeQpResidual() override;
/// UserObject that holds information (number of phases, components, etc)
const PorousFlowDictator & _dictator;
/// Number of fluid phases
const unsigned int _num_phases;
/// Whether fluid is present
const bool _fluid_present;
/// Porosity
const ADMaterialProperty<Real> & _porosity;
const MaterialProperty<Real> & _porosity_old;
/// Fluid density
const ADMaterialProperty<std::vector<Real>> * const _density;
const MaterialProperty<std::vector<Real>> * const _density_old;
/// Internal energy of porous matrix
const ADMaterialProperty<Real> & _rock_energy;
const MaterialProperty<Real> & _rock_energy_old;
/// Internal energy of fluid
const ADMaterialProperty<std::vector<Real>> * const _energy;
const MaterialProperty<std::vector<Real>> * const _energy_old;
/// Fluid phase saturation
const ADMaterialProperty<std::vector<Real>> * const _saturation;
const MaterialProperty<std::vector<Real>> * const _saturation_old;
};
(modules/navier_stokes/include/fvkernels/PCNSFVEnergyTimeDerivative.h)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#pragma once
#include "FVTimeKernel.h"
class PCNSFVEnergyTimeDerivative : public FVTimeKernel
{
public:
static InputParameters validParams();
PCNSFVEnergyTimeDerivative(const InputParameters & params);
protected:
ADReal computeQpResidual() override;
/// the density
const ADMaterialProperty<Real> & _rho;
/// the density time derivative
const ADMaterialProperty<Real> * _rho_dot;
/// the heat conductivity
const ADMaterialProperty<Real> & _cp;
/// the heat conductivity time derivative
const ADMaterialProperty<Real> * _cp_dot;
/// the porosity
const VariableValue & _eps;
/// whether this kernel is being used for a solid or a fluid temperature
const bool _is_solid;
/// scales the value of the kernel, used for faster steady state during pseudo transient
const Real _scaling;
/// whether a zero scaling factor has been specifed
const bool _zero_scaling;
};
(modules/navier_stokes/include/fvkernels/FVMatPropTimeKernel.h)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#pragma once
#include "FVTimeKernel.h"
/**
* Applies a residual equal to a supplied material property which is supposed to represent a time
* derivative, e.g. the derivative of the product of density and a passive scalar
* \f$\frac{\partial \rho s}{\partial t}\f$
*/
class FVMatPropTimeKernel : public FVTimeKernel
{
public:
static InputParameters validParams();
FVMatPropTimeKernel(const InputParameters & parameters);
protected:
ADReal computeQpResidual() override;
const ADMaterialProperty<Real> & _mat_prop_time_derivative;
};
(modules/heat_transfer/include/fvkernels/FVHeatConductionTimeDerivative.h)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#pragma once
#include "FVTimeKernel.h"
class FVHeatConductionTimeDerivative : public FVTimeKernel
{
public:
static InputParameters validParams();
FVHeatConductionTimeDerivative(const InputParameters & parameters);
protected:
virtual ADReal computeQpResidual() override;
/// Specific heat material property
const ADMaterialProperty<Real> & _specific_heat;
/// Density material property
const ADMaterialProperty<Real> & _density;
};