- group_fluxesAll the variables that hold the group fluxes. These MUST be listed by decreasing energy/increasing group number.
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:All the variables that hold the group fluxes. These MUST be listed by decreasing energy/increasing group number.
- num_groupsThe total numer of energy groups
C++ Type:unsigned int
Unit:(no unit assumed)
Controllable:No
Description:The total numer of energy groups
- powerThe reactor power.
C++ Type:double
Unit:(no unit assumed)
Controllable:No
Description:The reactor power.
- tot_fission_heatThe total fission heat postprocessor that's used to normalize the heat source.
C++ Type:PostprocessorName
Unit:(no unit assumed)
Controllable:No
Description:The total fission heat postprocessor that's used to normalize the heat source.
- 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
FissionHeatSource
Overview
This object adds the heat source term of the steady-state temperature advection-diffusion equation with coupled neutron fluxes from a -eigenvalue calculation. The power
input parameter allows the user to set a fixed thermal power output regardless of the neutron flux normalization factor. Alternatively, one may choose to normalize the neutron flux instead using the normalization
and normal_factor
parameters for the InversePowerMethod or NonlinearEigen executioners.
Example Input File Syntax
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
- displacementsThe displacements
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:The displacements
- 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_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.
- diag_save_inThe name of auxiliary variables to save this Kernel's diagonal Jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector<AuxVariableName>
Unit:(no unit assumed)
Controllable:No
Description:The name of auxiliary variables to save this Kernel's diagonal Jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
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
- save_inThe name of auxiliary variables to save this Kernel's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector<AuxVariableName>
Unit:(no unit assumed)
Controllable:No
Description:The name of auxiliary variables to save this Kernel's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
- seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
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
Input Files
- (problems/nts-temp-delayed.i)
- (problems/gmsh-two-mat-multiple-pin.i)
- (problems/gmsh-3d-coupled-nts-temp.i)
- (problems/msr-couple-temp-nts-2x2-10xPow.i)
- (problems/rodMatchSerpentWorth/rodded.i)
- (problems/rodMatchSerpentWorth/adjustAbsorb/rodded.i)
- (problems/gmsh-mesh-one-material.i)
- (tests/nts/gen-mesh-one-material.i)
- (problems/gmsh-mesh-two-materials.i)
- (problems/debug-precursor-action.i)
- (problems/axi-nts-temp-delayed.i)
- (problems/msr-couple-temp-nts-2x2.i)
- (problems/gen-mesh-one-material.i)
(problems/nts-temp-delayed.i)
[GlobalParams]
num_groups = 2
# temperature = temp
group_fluxes = 'group1 group2'
# MSRE full power = 10 MW; core volume 90 ft3
# power = 20
[]
[Mesh]
file = 'cylinder.msh'
[]
[Variables]
[group1]
order = FIRST
family = LAGRANGE
[]
[group2]
order = FIRST
family = LAGRANGE
[]
# [./temp]
# order = FIRST
# family = LAGRANGE
# # scaling = 1e-3
# [../]
[]
[Kernels]
# Neutronics
# [./time_group1]
# type = NtTimeDerivative
# grou_number = 1
# variable = group1
# [../]
# [./time_group2]
# type = NtTimeDerivative
# grou_number = 2
# variable = group2
# [../]
[diff_group1]
type = GroupDiffusion
variable = group1
group_number = 1
[]
[diff_group2]
type = GroupDiffusion
variable = group2
group_number = 2
[]
[sigma_r_group1]
type = SigmaR
variable = group1
group_number = 1
[]
[sigma_r_group2]
type = SigmaR
variable = group2
group_number = 2
[]
[inscatter_group1]
type = InScatter
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[]
[inscatter_group2]
type = InScatter
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[]
[fission_source_group1]
type = CoupledFissionEigenKernel
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[]
[fission_source_group2]
type = CoupledFissionEigenKernel
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[]
# [./fission_source_group1]
# type = CoupledFissionKernel
# variable = group1
# group_number = 1
# num_groups = 2
# group_fluxes = 'group1 group2'
# [../]
# [./fission_source_group2]
# type = CoupledFissionKernel
# variable = group2
# group_number = 2
# num_groups = 2
# group_fluxes = 'group1 group2'
# [../]
# # Temperature
# [./temp_cond]
# type = MatDiffusion
# variable = temp
# prop_name = 'k'
# save_in = 'diffus_resid tot_resid'
# [../]
# [./temp_source]
# type = FissionHeatSource
# tot_fissions = tot_fissions
# variable = temp
# save_in = 'src_resid tot_resid'
# [../]
[]
[AuxVariables]
# [./Qf]
# family = MONOMIAL
# order = CONSTANT
# [../]
# [./diffus_temp]
# family = MONOMIAL
# order = CONSTANT
# [../]
# [./diffus_resid]
# family = LAGRANGE
# order = FIRST
# [../]
# [./src_resid]
# family = LAGRANGE
# order = FIRST
# [../]
# [./bc_resid]
# family = LAGRANGE
# order = FIRST
# [../]
# [./tot_resid]
# family = LAGRANGE
# order = FIRST
# [../]
[]
[AuxKernels]
# [./Qf]
# type = FissionHeatSourceAux
# variable = Qf
# tot_fissions = tot_fissions
# [../]
# [./diffus_temp]
# type = MatDiffusionAux
# variable = diffus_temp
# diffuse_var = temp
# prop_name = 'k'
# [../]
[]
[Materials]
[fuel]
type = GenericMoltresMaterial
block = 'fuel'
property_tables_root = 'msr2g_enrU_mod_953_fuel_interp_'
num_groups = 2
num_precursor_groups = 8
prop_names = 'k'
prop_values = '.0123' # Cammi 2011 at 908 K
# prop_names = 'k d_k_d_temp'
# prop_values = '.0123 0' # Cammi 2011 at 908 K
[]
[moder]
type = GenericMoltresMaterial
block = 'moder'
property_tables_root = 'msr2g_enrU_fuel_922_mod_interp_'
num_groups = 2
num_precursor_groups = 8
prop_names = 'k'
prop_values = '.312' # Cammi 2011 at 908 K
# prop_names = 'k d_k_d_temp'
# prop_values = '.312 0' # Cammi 2011 at 908 K
[]
[]
[BCs]
# [./temp]
# boundary = boundary
# type = DirichletBC
# variable = temp
# value = 900
# save_in = 'bc_resid tot_resid'
# [../]
[group1_vacuum]
type = VacuumBC
variable = group1
boundary = 'all_top all_bottom'
[]
[group2_vacuum]
type = VacuumBC
variable = group2
boundary = 'all_top all_bottom'
[]
[]
[Executioner]
# type = NonlinearEigen
# free_power_iterations = 2
# source_abs_tol = 1e-12
# source_rel_tol = 1e-8
# output_after_power_iterations = false
type = InversePowerMethod
max_power_iterations = 50
xdiff = 'group1diff'
bx_norm = 'bnorm'
k0 = 2
pfactor = 1e-2
l_max_its = 100
# solve_type = 'PJFNK'
solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
petsc_options_iname = '-pc_type -sub_pc_type'
petsc_options_value = 'asm lu'
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Postprocessors]
[bnorm]
type = ElmIntegTotFissNtsPostprocessor
group_fluxes = 'group1 group2'
execute_on = linear
[]
[tot_fissions]
type = ElmIntegTotFissPostprocessor
execute_on = linear
[]
[group1norm]
type = ElementIntegralVariablePostprocessor
variable = group1
execute_on = linear
[]
[group2norm]
type = ElementIntegralVariablePostprocessor
variable = group2
execute_on = linear
[]
[group1max]
type = NodalExtremeValue
value_type = max
variable = group1
execute_on = timestep_end
[]
[group2max]
type = NodalExtremeValue
value_type = max
variable = group2
execute_on = timestep_end
[]
[group1diff]
type = ElementL2Diff
variable = group1
execute_on = 'linear timestep_end'
use_displaced_mesh = false
[]
[]
[Outputs]
[out]
type = Exodus
execute_on = 'timestep_end'
[]
[]
[Debug]
show_var_residual_norms = true
[]
# [ICs]
# [./temp_ic]
# type = ConstantIC
# variable = temp
# value = 900
# [../]
# []
(problems/gmsh-two-mat-multiple-pin.i)
[GlobalParams]
num_groups = 2
temperature = temp
group_fluxes = 'group1 group2'
# MSRE full power = 10 MW; core volume 90 ft3
power = 10
[]
[Mesh]
file = '/home/lindsayad/gdrive/gmsh-scripts/2d-msr-more-pins.msh'
[]
[Variables]
[group1]
order = FIRST
family = LAGRANGE
[]
[group2]
order = FIRST
family = LAGRANGE
[]
[temp]
order = FIRST
family = LAGRANGE
scaling = 1e-6
[]
[]
[Kernels]
# Neutronics
[diff_group1]
type = GroupDiffusion
variable = group1
group_number = 1
[]
[diff_group2]
type = GroupDiffusion
variable = group2
group_number = 2
[]
[sigma_r_group1]
type = SigmaR
variable = group1
group_number = 1
[]
[sigma_r_group2]
type = SigmaR
variable = group2
group_number = 2
[]
[inscatter_group1]
type = InScatter
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[]
[inscatter_group2]
type = InScatter
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[]
[fission_source_group1]
type = CoupledFissionEigenKernel
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[]
[fission_source_group2]
type = CoupledFissionEigenKernel
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[]
# Temperature
[temp_cond]
type = MatDiffusion
variable = temp
prop_name = 'k'
save_in = 'diffus_resid tot_resid'
[]
# [./temp_cond]
# type = HeatConduction
# diffusion_coefficient_name = k
# diffusion_coefficient_dT_name = d_k_d_temp
# use_displaced_mesh = false
# variable = temp
# [../]
[temp_source]
type = FissionHeatSource
tot_fissions = tot_fissions
variable = temp
save_in = 'src_resid tot_resid'
[]
[]
[AuxVariables]
[Qf]
family = MONOMIAL
order = CONSTANT
[]
[diffus_temp]
family = MONOMIAL
order = CONSTANT
[]
[diffus_resid]
family = LAGRANGE
order = FIRST
[]
[src_resid]
family = LAGRANGE
order = FIRST
[]
[bc_resid]
family = LAGRANGE
order = FIRST
[]
[tot_resid]
family = LAGRANGE
order = FIRST
[]
[]
[AuxKernels]
[Qf]
type = FissionHeatSourceAux
variable = Qf
tot_fissions = tot_fissions
[]
[diffus_temp]
type = MatDiffusionAux
variable = diffus_temp
diffuse_var = temp
prop_name = 'k'
[]
[]
[Materials]
[fuel]
type = GenericMoltresMaterial
block = 'fuel'
property_tables_root = '/home/lindsayad/serpent/core/examples/serpent-input/msre/msr2g_enrU_mod_953_fuel_interp_'
num_groups = 2
prop_names = 'k'
prop_values = '.0123' # Cammi 2011 at 908 K
# prop_names = 'k d_k_d_temp'
# prop_values = '.0123 0' # Cammi 2011 at 908 K
[]
[moder]
type = GenericMoltresMaterial
block = 'moder'
property_tables_root = '/home/lindsayad/serpent/core/examples/serpent-input/msre/msr2g_enrU_fuel_922_mod_interp_'
num_groups = 2
prop_names = 'k'
prop_values = '.312' # Cammi 2011 at 908 K
# prop_names = 'k d_k_d_temp'
# prop_values = '.312 0' # Cammi 2011 at 908 K
[]
[]
[BCs]
[temp]
boundary = 'boundary'
type = DirichletBC
variable = temp
value = 900
save_in = 'bc_resid tot_resid'
[]
# [./temp]
# boundary = boundary
# type = VacuumBC
# variable = temp
# [../]
[]
[Executioner]
# type = NonlinearEigen
# free_power_iterations = 4
# source_abs_tol = 1e-12
# source_rel_tol = 1e-8
# output_after_power_iterations = true
type = InversePowerMethod
max_power_iterations = 50
xdiff = 'group1diff'
bx_norm = 'bnorm'
k0 = 1.0
pfactor = 1e-2
l_max_its = 100
solve_type = 'PJFNK'
# solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
# petsc_options_iname = '-pc_type -sub_pc_type'
# petsc_options_value = 'asm lu'
[]
[Preconditioning]
[SMP]
type = SMP
# full = true
[]
[]
[Postprocessors]
[bnorm]
type = ElmIntegTotFissNtsPostprocessor
group_fluxes = 'group1 group2'
execute_on = linear
[]
[tot_fissions]
type = ElmIntegTotFissPostprocessor
execute_on = linear
[]
[group1norm]
type = ElementIntegralVariablePostprocessor
variable = group1
execute_on = linear
[]
[group2norm]
type = ElementIntegralVariablePostprocessor
variable = group2
execute_on = linear
[]
[group1max]
type = NodalExtremeValue
value_type = max
variable = group1
execute_on = timestep_end
[]
[group2max]
type = NodalExtremeValue
value_type = max
variable = group2
execute_on = timestep_end
[]
[src_resid_post]
type = NodalL2Norm
variable = src_resid
execute_on = nonlinear
[]
[diffus_resid_post]
type = NodalL2Norm
variable = diffus_resid
execute_on = nonlinear
[]
[bc_resid_post]
type = NodalL2Norm
variable = bc_resid
execute_on = nonlinear
[]
[tot_resid_post]
type = NodalL2Norm
variable = tot_resid
execute_on = nonlinear
[]
[group1diff]
type = ElementL2Diff
variable = group1
execute_on = 'linear timestep_end'
use_displaced_mesh = false
[]
[]
[Outputs]
[out]
type = Exodus
execute_on = 'timestep_end'
[]
[]
[Debug]
show_var_residual_norms = true
[]
[ICs]
[temp_ic]
type = ConstantIC
variable = temp
value = 900
[]
[]
(problems/gmsh-3d-coupled-nts-temp.i)
[GlobalParams]
num_groups = 2
temperature = temp
group_fluxes = 'group1 group2'
# MSRE full power = 10 MW; core volume 90 ft3
power = 20
[]
[Mesh]
file = 'msr-small.msh'
[]
[Variables]
[group1]
order = FIRST
family = LAGRANGE
[]
[group2]
order = FIRST
family = LAGRANGE
[]
[temp]
order = FIRST
family = LAGRANGE
# scaling = 1e-3
[]
[]
[Kernels]
# Neutronics
[diff_group1]
type = GroupDiffusion
variable = group1
group_number = 1
[]
[diff_group2]
type = GroupDiffusion
variable = group2
group_number = 2
[]
[sigma_r_group1]
type = SigmaR
variable = group1
group_number = 1
[]
[sigma_r_group2]
type = SigmaR
variable = group2
group_number = 2
[]
[inscatter_group1]
type = InScatter
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[]
[inscatter_group2]
type = InScatter
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[]
[fission_source_group1]
type = CoupledFissionEigenKernel
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[]
[fission_source_group2]
type = CoupledFissionEigenKernel
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[]
# Temperature
[temp_cond]
type = MatDiffusion
variable = temp
prop_name = 'k'
save_in = 'diffus_resid tot_resid'
[]
# [./temp_cond]
# type = HeatConduction
# diffusion_coefficient_name = k
# diffusion_coefficient_dT_name = d_k_d_temp
# use_displaced_mesh = false
# variable = temp
# [../]
[temp_source]
type = FissionHeatSource
tot_fissions = tot_fissions
variable = temp
save_in = 'src_resid tot_resid'
[]
[]
[AuxVariables]
[Qf]
family = MONOMIAL
order = CONSTANT
[]
[diffus_temp]
family = MONOMIAL
order = CONSTANT
[]
[diffus_resid]
family = LAGRANGE
order = FIRST
[]
[src_resid]
family = LAGRANGE
order = FIRST
[]
[bc_resid]
family = LAGRANGE
order = FIRST
[]
[tot_resid]
family = LAGRANGE
order = FIRST
[]
[]
[AuxKernels]
[Qf]
type = FissionHeatSourceAux
variable = Qf
tot_fissions = tot_fissions
[]
[diffus_temp]
type = MatDiffusionAux
variable = diffus_temp
diffuse_var = temp
prop_name = 'k'
[]
[]
[Materials]
[fuel]
type = GenericMoltresMaterial
block = 'fuel'
property_tables_root = 'msr2g_enrU_mod_953_fuel_interp_'
num_groups = 2
prop_names = 'k'
prop_values = '.0123' # Cammi 2011 at 908 K
# prop_names = 'k d_k_d_temp'
# prop_values = '.0123 0' # Cammi 2011 at 908 K
[]
[moder]
type = GenericMoltresMaterial
block = 'moder'
property_tables_root = 'msr2g_enrU_fuel_922_mod_interp_'
num_groups = 2
prop_names = 'k'
prop_values = '.312' # Cammi 2011 at 908 K
# prop_names = 'k d_k_d_temp'
# prop_values = '.312 0' # Cammi 2011 at 908 K
[]
[]
[BCs]
[temp]
boundary = boundary
type = DirichletBC
variable = temp
value = 900
save_in = 'bc_resid tot_resid'
[]
# [./temp]
# boundary = boundary
# type = VacuumBC
# variable = temp
# [../]
[]
[Executioner]
type = NonlinearEigen
free_power_iterations = 2
source_abs_tol = 1e-12
source_rel_tol = 1e-8
output_after_power_iterations = false
# type = InversePowerMethod
# max_power_iterations = 50
# xdiff = 'group1diff'
bx_norm = 'bnorm'
k0 = 2
pfactor = 1e-2
l_max_its = 100
solve_type = 'PJFNK'
# solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
petsc_options_iname = '-pc_type -sub_pc_type'
petsc_options_value = 'asm lu'
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Postprocessors]
[bnorm]
type = ElmIntegTotFissNtsPostprocessor
group_fluxes = 'group1 group2'
execute_on = linear
[]
[tot_fissions]
type = ElmIntegTotFissPostprocessor
execute_on = linear
[]
[group1norm]
type = ElementIntegralVariablePostprocessor
variable = group1
execute_on = linear
[]
[group2norm]
type = ElementIntegralVariablePostprocessor
variable = group2
execute_on = linear
[]
[group1max]
type = NodalExtremeValue
value_type = max
variable = group1
execute_on = timestep_end
[]
[group2max]
type = NodalExtremeValue
value_type = max
variable = group2
execute_on = timestep_end
[]
# [./src_resid_post]
# type = NodalL2Norm
# variable = src_resid
# execute_on = nonlinear
# [../]
# [./diffus_resid_post]
# type = NodalL2Norm
# variable = diffus_resid
# execute_on = nonlinear
# [../]
# [./bc_resid_post]
# type = NodalL2Norm
# variable = bc_resid
# execute_on = nonlinear
# [../]
# [./tot_resid_post]
# type = NodalL2Norm
# variable = tot_resid
# execute_on = nonlinear
# [../]
# [./group1diff]
# type = ElementL2Diff
# variable = group1
# execute_on = 'linear timestep_end'
# use_displaced_mesh = false
# [../]
[]
[Outputs]
[out]
type = Exodus
execute_on = 'timestep_end'
[]
[]
[Debug]
show_var_residual_norms = true
[]
# [ICs]
# [./temp_ic]
# type = ConstantIC
# variable = temp
# value = 900
# [../]
# []
(problems/msr-couple-temp-nts-2x2-10xPow.i)
[GlobalParams]
num_groups = 2
temperature = temp
group_fluxes = 'group1 group2'
# MSRE full power = 10 MW; core volume 90 ft3
power = 20
[]
[Mesh]
file = '/home/lindsayad/gdrive/gmsh-scripts/msr-small.msh'
[]
[Variables]
[group1]
order = FIRST
family = LAGRANGE
[]
[group2]
order = FIRST
family = LAGRANGE
[]
[temp]
order = FIRST
family = LAGRANGE
# scaling = 1e-3
[]
[]
[Kernels]
# Neutronics
[diff_group1]
type = GroupDiffusion
variable = group1
group_number = 1
[]
[diff_group2]
type = GroupDiffusion
variable = group2
group_number = 2
[]
[sigma_r_group1]
type = SigmaR
variable = group1
group_number = 1
[]
[sigma_r_group2]
type = SigmaR
variable = group2
group_number = 2
[]
[inscatter_group1]
type = InScatter
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[]
[inscatter_group2]
type = InScatter
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[]
[fission_source_group1]
type = CoupledFissionEigenKernel
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[]
[fission_source_group2]
type = CoupledFissionEigenKernel
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[]
# Temperature
[temp_cond]
type = MatDiffusion
variable = temp
prop_name = 'k'
save_in = 'diffus_resid tot_resid'
[]
# [./temp_cond]
# type = HeatConduction
# diffusion_coefficient_name = k
# diffusion_coefficient_dT_name = d_k_d_temp
# use_displaced_mesh = false
# variable = temp
# [../]
[temp_source]
type = FissionHeatSource
tot_fissions = tot_fissions
variable = temp
save_in = 'src_resid tot_resid'
[]
[]
[AuxVariables]
[Qf]
family = MONOMIAL
order = CONSTANT
[]
[diffus_temp]
family = MONOMIAL
order = CONSTANT
[]
[diffus_resid]
family = LAGRANGE
order = FIRST
[]
[src_resid]
family = LAGRANGE
order = FIRST
[]
[bc_resid]
family = LAGRANGE
order = FIRST
[]
[tot_resid]
family = LAGRANGE
order = FIRST
[]
[]
[AuxKernels]
[Qf]
type = FissionHeatSourceAux
variable = Qf
tot_fissions = tot_fissions
[]
[diffus_temp]
type = MatDiffusionAux
variable = diffus_temp
diffuse_var = temp
prop_name = 'k'
[]
[]
[Materials]
[fuel]
type = GenericMoltresMaterial
block = 'fuel'
property_tables_root = '/home/lindsayad/serpent/core/examples/serpent-input/msre/msr2g_enrU_mod_953_fuel_interp_'
num_groups = 2
prop_names = 'k'
prop_values = '.0123' # Cammi 2011 at 908 K
# prop_names = 'k d_k_d_temp'
# prop_values = '.0123 0' # Cammi 2011 at 908 K
[]
[moder]
type = GenericMoltresMaterial
block = 'moder'
property_tables_root = '/home/lindsayad/serpent/core/examples/serpent-input/msre/msr2g_enrU_fuel_922_mod_interp_'
num_groups = 2
prop_names = 'k'
prop_values = '.312' # Cammi 2011 at 908 K
# prop_names = 'k d_k_d_temp'
# prop_values = '.312 0' # Cammi 2011 at 908 K
[]
[]
[BCs]
[temp]
boundary = boundary
type = DirichletBC
variable = temp
value = 900
save_in = 'bc_resid tot_resid'
[]
# [./temp]
# boundary = boundary
# type = VacuumBC
# variable = temp
# [../]
[]
[Executioner]
type = NonlinearEigen
free_power_iterations = 2
source_abs_tol = 1e-12
source_rel_tol = 1e-8
output_after_power_iterations = false
# type = InversePowerMethod
# max_power_iterations = 50
# xdiff = 'group1diff'
bx_norm = 'bnorm'
k0 = 1.0
pfactor = 1e-2
l_max_its = 100
# solve_type = 'PJFNK'
solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
petsc_options_iname = '-pc_type -sub_pc_type'
petsc_options_value = 'asm lu'
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Postprocessors]
[bnorm]
type = ElmIntegTotFissNtsPostprocessor
group_fluxes = 'group1 group2'
execute_on = linear
[]
[tot_fissions]
type = ElmIntegTotFissPostprocessor
execute_on = linear
[]
[group1norm]
type = ElementIntegralVariablePostprocessor
variable = group1
execute_on = linear
[]
[group2norm]
type = ElementIntegralVariablePostprocessor
variable = group2
execute_on = linear
[]
[group1max]
type = NodalExtremeValue
value_type = max
variable = group1
execute_on = timestep_end
[]
[group2max]
type = NodalExtremeValue
value_type = max
variable = group2
execute_on = timestep_end
[]
[src_resid_post]
type = NodalL2Norm
variable = src_resid
execute_on = nonlinear
[]
[diffus_resid_post]
type = NodalL2Norm
variable = diffus_resid
execute_on = nonlinear
[]
[bc_resid_post]
type = NodalL2Norm
variable = bc_resid
execute_on = nonlinear
[]
[tot_resid_post]
type = NodalL2Norm
variable = tot_resid
execute_on = nonlinear
[]
[group1diff]
type = ElementL2Diff
variable = group1
execute_on = 'linear timestep_end'
use_displaced_mesh = false
[]
[]
[Outputs]
[out]
type = Exodus
execute_on = 'timestep_end'
[]
[]
[Debug]
show_var_residual_norms = true
[]
# [ICs]
# [./temp_ic]
# type = ConstantIC
# variable = temp
# value = 900
# [../]
# []
(problems/rodMatchSerpentWorth/rodded.i)
flow_velocity = 21.7 # cm/s. See MSRE-properties.ods
nt_scale = 1e13
ini_temp = 922
diri_temp = 922
[GlobalParams]
num_groups = 4
num_precursor_groups = 6
use_exp_form = false
group_fluxes = 'group1 group2 group3 group4'
temperature = temp
sss2_input = true
pre_concs = 'pre1 pre2 pre3 pre4 pre5 pre6'
account_delayed = false
[]
[Mesh]
file = '2d_rodded_lattice.msh'
[]
[Problem]
coord_type = RZ
kernel_coverage_check = false
[]
[Variables]
[temp]
initial_condition = ${ini_temp}
scaling = 1e-4
[]
[]
[Precursors]
[pres]
var_name_base = pre
block = 'fuel'
outlet_boundaries = 'fuel_tops'
u_def = 0
v_def = ${flow_velocity}
w_def = 0
nt_exp_form = false
family = MONOMIAL
order = CONSTANT
# jac_test = true
[]
[]
[Nt]
var_name_base = group
vacuum_boundaries = 'fuel_bottoms fuel_tops moder_bottoms moder_tops outer_wall cRod_top cRod_bot'
create_temperature_var = false
eigen = true
[]
[Kernels]
[temp_source_fuel]
type = FissionHeatSource
variable = temp
nt_scale = ${nt_scale}
block = 'fuel'
power = 7.5e6
tot_fissions = tot_fissions
[]
[temp_source_mod]
type = GammaHeatSource
variable = temp
gamma = .0144 # Cammi .0144
block = 'moder'
average_fission_heat = 'tot_fissions'
[]
[temp_diffusion]
type = MatDiffusion
diffusivity = 'k'
variable = temp
[]
[temp_advection_fuel]
type = ConservativeTemperatureAdvection
velocity = '0 ${flow_velocity} 0'
variable = temp
block = 'fuel'
[]
[]
[BCs]
[temp_diri_cg]
boundary = 'moder_bottoms fuel_bottoms outer_wall'
type = DirichletBC
variable = temp
value = ${diri_temp}
[]
[temp_advection_outlet]
boundary = 'fuel_tops'
type = TemperatureOutflowBC
variable = temp
velocity = '0 ${flow_velocity} 0'
[]
[]
[Materials]
[fuel]
type = GenericMoltresMaterial
property_tables_root = '../../tutorial/step01_groupConstants/MSREProperties/msre_gentry_4gfuel_'
interp_type = 'spline'
block = 'fuel'
prop_names = 'k cp'
prop_values = '.0553 1967' # Robertson MSRE technical report @ 922 K
[]
[rho_fuel]
type = DerivativeParsedMaterial
f_name = rho
function = '2.146e-3 * exp(-1.8 * 1.18e-4 * (temp - 922))'
args = 'temp'
derivative_order = 1
block = 'fuel'
[]
[moder]
type = GenericMoltresMaterial
property_tables_root = '../../tutorial/step01_groupConstants/MSREProperties/msre_gentry_4gmoder_'
interp_type = 'spline'
prop_names = 'k cp'
prop_values = '.312 1760' # Cammi 2011 at 908 K
block = 'moder'
[]
[rho_moder]
type = DerivativeParsedMaterial
f_name = rho
function = '1.86e-3 * exp(-1.8 * 1.0e-5 * (temp - 922))'
args = 'temp'
derivative_order = 1
block = 'moder'
[]
[cRod]
type = RoddedMaterial
property_tables_root = '../../tutorial/step01_groupConstants/MSREProperties/msre_gentry_4gmoder_'
interp_type = 'spline'
prop_names = 'k cp'
prop_values = '.312 1760' # Cammi 2011 at 908 K
block = 'cRod'
rodDimension = 'y'
rodPosition = 200.0
absorb_factor = 1e6 # how much more absorbing than usual in absorbing region?
[]
[rho_crod]
type = DerivativeParsedMaterial
f_name = rho
function = '1.86e-3 * exp(-1.8 * 1.0e-5 * (temp - 922))'
args = 'temp'
derivative_order = 1
block = 'cRod'
[]
[]
[Executioner]
type = InversePowerMethod
max_power_iterations = 50
xdiff = 'group1diff'
bx_norm = 'bnorm'
k0 = 1.5
pfactor = 1e-2
l_max_its = 100
# solve_type = 'PJFNK'
solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
petsc_options_iname = '-pc_type -sub_pc_type -pc_factor_shift_type -pc_factor_shift_amount'
petsc_options_value = 'asm lu NONZERO 1e-10'
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Postprocessors]
[bnorm]
type = ElmIntegTotFissNtsPostprocessor
execute_on = linear
[]
[tot_fissions]
type = ElmIntegTotFissPostprocessor
execute_on = linear
[]
[group1norm]
type = ElementIntegralVariablePostprocessor
variable = group1
execute_on = linear
[]
[group1max]
type = NodalExtremeValue
value_type = max
variable = group1
execute_on = timestep_end
[]
[group1diff]
type = ElementL2Diff
variable = group1
execute_on = 'linear timestep_end'
use_displaced_mesh = false
[]
[group2norm]
type = ElementIntegralVariablePostprocessor
variable = group2
execute_on = linear
[]
[group2max]
type = NodalExtremeValue
value_type = max
variable = group2
execute_on = timestep_end
[]
[group2diff]
type = ElementL2Diff
variable = group2
execute_on = 'linear timestep_end'
use_displaced_mesh = false
[]
[]
[Outputs]
perf_graph = true
print_linear_residuals = true
csv = true
exodus = true
execute_on = 'final'
file_base = 'out'
[]
[Debug]
show_var_residual_norms = true
[]
(problems/rodMatchSerpentWorth/adjustAbsorb/rodded.i)
flow_velocity = 21.7 # cm/s. See MSRE-properties.ods
nt_scale = 1e13
ini_temp = 922
diri_temp = 922
[GlobalParams]
num_groups = 4
num_precursor_groups = 6
use_exp_form = false
group_fluxes = 'group1 group2 group3 group4'
temperature = temp
sss2_input = true
pre_concs = 'pre1 pre2 pre3 pre4 pre5 pre6'
account_delayed = false
[]
[Mesh]
file = '2d_rodded_lattice.msh'
[]
[Problem]
coord_type = RZ
kernel_coverage_check = false
[]
[Variables]
[temp]
initial_condition = ${ini_temp}
scaling = 1e-4
[]
[]
[Precursors]
[pres]
var_name_base = pre
block = 'fuel'
outlet_boundaries = 'fuel_tops'
u_def = 0
v_def = ${flow_velocity}
w_def = 0
nt_exp_form = false
family = MONOMIAL
order = CONSTANT
# jac_test = true
[]
[]
[Nt]
var_name_base = group
vacuum_boundaries = 'fuel_bottoms fuel_tops moder_bottoms moder_tops outer_wall cRod_top cRod_bot'
create_temperature_var = false
eigen = true
[]
[Kernels]
[temp_source_fuel]
type = FissionHeatSource
variable = temp
nt_scale = ${nt_scale}
block = 'fuel'
power = 7.5e6
tot_fissions = tot_fissions
[]
[temp_source_mod]
type = GammaHeatSource
variable = temp
gamma = .0144 # Cammi .0144
block = 'moder'
average_fission_heat = 'tot_fissions'
[]
[temp_diffusion]
type = MatDiffusion
diffusivity = 'k'
variable = temp
[]
[temp_advection_fuel]
type = ConservativeTemperatureAdvection
velocity = '0 ${flow_velocity} 0'
variable = temp
block = 'fuel'
[]
[]
[BCs]
[temp_diri_cg]
boundary = 'moder_bottoms fuel_bottoms outer_wall'
type = DirichletBC
variable = temp
value = ${diri_temp}
[]
[temp_advection_outlet]
boundary = 'fuel_tops'
type = TemperatureOutflowBC
variable = temp
velocity = '0 ${flow_velocity} 0'
[]
[]
[Materials]
[fuel]
type = GenericMoltresMaterial
property_tables_root = '../../../tutorial/step01_groupConstants/MSREProperties/msre_gentry_4gfuel_'
interp_type = 'spline'
block = 'fuel'
prop_names = 'k cp'
prop_values = '.0553 1967' # Robertson MSRE technical report @ 922 K
[]
[rho_fuel]
type = DerivativeParsedMaterial
f_name = rho
function = '2.146e-3 * exp(-1.8 * 1.18e-4 * (temp - 922))'
args = 'temp'
derivative_order = 1
block = 'fuel'
[]
[moder]
type = GenericMoltresMaterial
property_tables_root = '../../../tutorial/step01_groupConstants/MSREProperties/msre_gentry_4gmoder_'
interp_type = 'spline'
prop_names = 'k cp'
prop_values = '.312 1760' # Cammi 2011 at 908 K
block = 'moder'
[]
[rho_moder]
type = DerivativeParsedMaterial
f_name = rho
function = '1.86e-3 * exp(-1.8 * 1.0e-5 * (temp - 922))'
args = 'temp'
derivative_order = 1
block = 'moder'
[]
[cRod]
type = RoddedMaterial
property_tables_root = '../../../tutorial/step01_groupConstants/MSREProperties/msre_gentry_4gmoder_'
interp_type = 'spline'
prop_names = 'k cp'
prop_values = '.312 1760' # Cammi 2011 at 908 K
block = 'cRod'
rodDimension = 'y'
rodPosition = 23.622
absorb_factor = 14.22 # how much more absorbing than usual in absorbing region?
[]
[rho_crod]
type = DerivativeParsedMaterial
f_name = rho
function = '1.86e-3 * exp(-1.8 * 1.0e-5 * (temp - 922))'
args = 'temp'
derivative_order = 1
block = 'cRod'
[]
[]
[Executioner]
type = InversePowerMethod
max_power_iterations = 50
xdiff = 'group1diff'
bx_norm = 'bnorm'
k0 = 1.5
pfactor = 1e-2
l_max_its = 100
# solve_type = 'PJFNK'
solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
petsc_options_iname = '-pc_type -sub_pc_type -pc_factor_shift_type -pc_factor_shift_amount'
petsc_options_value = 'asm lu NONZERO 1e-10'
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Postprocessors]
[bnorm]
type = ElmIntegTotFissNtsPostprocessor
execute_on = linear
[]
[tot_fissions]
type = ElmIntegTotFissPostprocessor
execute_on = linear
[]
[group1norm]
type = ElementIntegralVariablePostprocessor
variable = group1
execute_on = linear
[]
[group1max]
type = NodalExtremeValue
value_type = max
variable = group1
execute_on = timestep_end
[]
[group1diff]
type = ElementL2Diff
variable = group1
execute_on = 'linear timestep_end'
use_displaced_mesh = false
[]
[group2norm]
type = ElementIntegralVariablePostprocessor
variable = group2
execute_on = linear
[]
[group2max]
type = NodalExtremeValue
value_type = max
variable = group2
execute_on = timestep_end
[]
[group2diff]
type = ElementL2Diff
variable = group2
execute_on = 'linear timestep_end'
use_displaced_mesh = false
[]
[]
[Outputs]
perf_graph = true
print_linear_residuals = true
csv = true
exodus = true
execute_on = 'final'
file_base = 'out'
[]
[Debug]
show_var_residual_norms = true
[]
(problems/gmsh-mesh-one-material.i)
[GlobalParams]
num_groups = 2
temperature = temp
group_fluxes = 'group1 group2'
# MSRE full power = 10 MW; core volume 90 ft3
power = 10
[]
[Mesh]
file = '/home/lindsayad/gdrive/gmsh-scripts/2d-one-mat.msh'
[]
[Variables]
[group1]
order = FIRST
family = LAGRANGE
[]
[group2]
order = FIRST
family = LAGRANGE
[]
[temp]
order = FIRST
family = LAGRANGE
scaling = 1e-6
[]
[]
[Kernels]
# Neutronics
[diff_group1]
type = GroupDiffusion
variable = group1
group_number = 1
[]
[diff_group2]
type = GroupDiffusion
variable = group2
group_number = 2
[]
[sigma_r_group1]
type = SigmaR
variable = group1
group_number = 1
[]
[sigma_r_group2]
type = SigmaR
variable = group2
group_number = 2
[]
[inscatter_group1]
type = InScatter
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[]
[inscatter_group2]
type = InScatter
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[]
[fission_source_group1]
type = CoupledFissionEigenKernel
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[]
[fission_source_group2]
type = CoupledFissionEigenKernel
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[]
# Temperature
[temp_cond]
type = MatDiffusion
variable = temp
prop_name = 'k'
save_in = 'diffus_resid tot_resid'
[]
# [./temp_cond]
# type = HeatConduction
# diffusion_coefficient_name = k
# diffusion_coefficient_dT_name = d_k_d_temp
# use_displaced_mesh = false
# variable = temp
# [../]
[temp_source]
type = FissionHeatSource
tot_fissions = tot_fissions
variable = temp
save_in = 'src_resid tot_resid'
[]
[]
[AuxVariables]
[Qf]
family = MONOMIAL
order = CONSTANT
[]
[diffus_temp]
family = MONOMIAL
order = CONSTANT
[]
[diffus_resid]
family = LAGRANGE
order = FIRST
[]
[src_resid]
family = LAGRANGE
order = FIRST
[]
[bc_resid]
family = LAGRANGE
order = FIRST
[]
[tot_resid]
family = LAGRANGE
order = FIRST
[]
[]
[AuxKernels]
[Qf]
type = FissionHeatSourceAux
variable = Qf
tot_fissions = tot_fissions
[]
[diffus_temp]
type = MatDiffusionAux
variable = diffus_temp
diffuse_var = temp
prop_name = 'k'
[]
[]
[Materials]
[fuel]
type = GenericMoltresMaterial
block = 0
property_tables_root = '/home/lindsayad/serpent/core/examples/serpent-input/msre/msr2g_enrU_mod_953_fuel_interp_'
num_groups = 2
prop_names = 'k'
prop_values = '.0123' # Cammi 2011 at 908 K
# prop_names = 'k d_k_d_temp'
# prop_values = '.0123 0' # Cammi 2011 at 908 K
[]
[]
[BCs]
[temp]
boundary = 'boundary'
type = DirichletBC
variable = temp
value = 900
save_in = 'bc_resid tot_resid'
[]
# [./temp]
# boundary = boundary
# type = VacuumBC
# variable = temp
# [../]
[]
[Executioner]
# type = NonlinearEigen
# free_power_iterations = 4
# source_abs_tol = 1e-12
# source_rel_tol = 1e-8
# output_after_power_iterations = true
type = InversePowerMethod
max_power_iterations = 50
xdiff = 'group1diff'
bx_norm = 'bnorm'
k0 = 1.0
pfactor = 1e-2
l_max_its = 100
solve_type = 'PJFNK'
# solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
# petsc_options_iname = '-pc_type -sub_pc_type'
# petsc_options_value = 'asm lu'
[]
[Preconditioning]
[SMP]
type = SMP
# full = true
[]
[]
[Postprocessors]
[bnorm]
type = ElmIntegTotFissNtsPostprocessor
group_fluxes = 'group1 group2'
execute_on = linear
[]
[tot_fissions]
type = ElmIntegTotFissPostprocessor
execute_on = linear
[]
[group1norm]
type = ElementIntegralVariablePostprocessor
variable = group1
execute_on = linear
[]
[group2norm]
type = ElementIntegralVariablePostprocessor
variable = group2
execute_on = linear
[]
[group1max]
type = NodalExtremeValue
value_type = max
variable = group1
execute_on = timestep_end
[]
[group2max]
type = NodalExtremeValue
value_type = max
variable = group2
execute_on = timestep_end
[]
[src_resid_post]
type = NodalL2Norm
variable = src_resid
execute_on = nonlinear
[]
[diffus_resid_post]
type = NodalL2Norm
variable = diffus_resid
execute_on = nonlinear
[]
[bc_resid_post]
type = NodalL2Norm
variable = bc_resid
execute_on = nonlinear
[]
[tot_resid_post]
type = NodalL2Norm
variable = tot_resid
execute_on = nonlinear
[]
[group1diff]
type = ElementL2Diff
variable = group1
execute_on = 'linear timestep_end'
use_displaced_mesh = false
[]
[]
[Outputs]
[out]
type = Exodus
execute_on = 'timestep_end'
[]
[]
[Debug]
show_var_residual_norms = true
[]
[ICs]
[temp_ic]
type = ConstantIC
variable = temp
value = 900
[]
[]
(tests/nts/gen-mesh-one-material.i)
[GlobalParams]
num_groups = 2
num_precursor_groups = 8
temperature = temp
group_fluxes = 'group1 group2'
# MSRE full power = 10 MW; core volume 90 ft3
power = 10
use_exp_form = false
sss2_input = false
account_delayed = false
[]
[Mesh]
type = GeneratedMesh
dim = 2
xmax = 6
ymax = 6
nx = 15
ny = 15
[]
[Variables]
[group1]
order = FIRST
family = LAGRANGE
[]
[group2]
order = FIRST
family = LAGRANGE
[]
[temp]
order = FIRST
family = LAGRANGE
scaling = 1e-6
[]
[]
[Kernels]
# Neutronics
[diff_group1]
type = GroupDiffusion
variable = group1
group_number = 1
[]
[diff_group2]
type = GroupDiffusion
variable = group2
group_number = 2
[]
[sigma_r_group1]
type = SigmaR
variable = group1
group_number = 1
[]
[sigma_r_group2]
type = SigmaR
variable = group2
group_number = 2
[]
[inscatter_group1]
type = InScatter
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[]
[inscatter_group2]
type = InScatter
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[]
[fission_source_group1]
type = CoupledFissionEigenKernel
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[]
[fission_source_group2]
type = CoupledFissionEigenKernel
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[]
# Temperature
[temp_cond]
type = MatDiffusion
variable = temp
diffusivity = 'k'
save_in = 'diffus_resid tot_resid'
[]
# [./temp_cond]
# type = HeatConduction
# diffusion_coefficient_name = k
# diffusion_coefficient_dT_name = d_k_d_temp
# use_displaced_mesh = false
# variable = temp
# [../]
[temp_source]
type = FissionHeatSource
tot_fission_heat = tot_fission_heat
variable = temp
save_in = 'src_resid tot_resid'
[]
[]
[AuxVariables]
[Qf]
family = MONOMIAL
order = CONSTANT
[]
[diffus_temp]
family = MONOMIAL
order = CONSTANT
[]
[diffus_resid]
family = LAGRANGE
order = FIRST
[]
[src_resid]
family = LAGRANGE
order = FIRST
[]
[bc_resid]
family = LAGRANGE
order = FIRST
[]
[tot_resid]
family = LAGRANGE
order = FIRST
[]
[]
[AuxKernels]
[Qf]
type = FissionHeatSourceAux
variable = Qf
tot_fission_heat = tot_fission_heat
[]
[diffus_temp]
type = MatDiffusionAux
variable = diffus_temp
diffuse_var = temp
prop_name = 'k'
[]
[]
[Materials]
[fuel]
type = GenericMoltresMaterial
block = 0
property_tables_root = 'msr2g_enrU_mod_953_fuel_interp_'
prop_names = 'k'
prop_values = '.0123' # Cammi 2011 at 908 K
interp_type = spline
[]
[]
[BCs]
[temp]
boundary = 'left right top bottom'
type = DirichletBC
variable = temp
value = 900
save_in = 'bc_resid tot_resid'
[]
# [./temp]
# boundary = boundary
# type = VacuumBC
# variable = temp
# [../]
[]
[Executioner]
# type = NonlinearEigen
# free_power_iterations = 4
# source_abs_tol = 1e-12
# source_rel_tol = 1e-8
# output_after_power_iterations = true
type = InversePowerMethod
max_power_iterations = 50
xdiff = 'group1diff'
bx_norm = 'bnorm'
k0 = 1.0
pfactor = 1e-2
l_max_its = 100
solve_type = 'PJFNK'
# solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
# petsc_options_iname = '-pc_type -sub_pc_type'
# petsc_options_value = 'asm lu'
[]
[Preconditioning]
[SMP]
type = SMP
# full = true
[]
[]
[Postprocessors]
[bnorm]
type = ElmIntegTotFissNtsPostprocessor
group_fluxes = 'group1 group2'
execute_on = linear
[]
[tot_fission_heat]
type = ElmIntegTotFissHeatPostprocessor
execute_on = linear
[]
[group1norm]
type = ElementIntegralVariablePostprocessor
variable = group1
execute_on = linear
[]
[group2norm]
type = ElementIntegralVariablePostprocessor
variable = group2
execute_on = linear
[]
[group1max]
type = NodalExtremeValue
value_type = max
variable = group1
execute_on = timestep_end
[]
[group2max]
type = NodalExtremeValue
value_type = max
variable = group2
execute_on = timestep_end
[]
[src_resid_post]
type = NodalL2Norm
variable = src_resid
execute_on = nonlinear
[]
[diffus_resid_post]
type = NodalL2Norm
variable = diffus_resid
execute_on = nonlinear
[]
[bc_resid_post]
type = NodalL2Norm
variable = bc_resid
execute_on = nonlinear
[]
[tot_resid_post]
type = NodalL2Norm
variable = tot_resid
execute_on = nonlinear
[]
[group1diff]
type = ElementL2Diff
variable = group1
execute_on = 'linear timestep_end'
use_displaced_mesh = false
[]
[]
[Outputs]
[out]
type = Exodus
execute_on = 'timestep_end'
[]
[]
[Debug]
show_var_residual_norms = true
[]
[ICs]
[temp_ic]
type = ConstantIC
variable = temp
value = 900
[]
[]
(problems/gmsh-mesh-two-materials.i)
[GlobalParams]
num_groups = 2
temperature = temp
group_fluxes = 'group1 group2'
# MSRE full power = 10 MW; core volume 90 ft3
power = 10
[]
[Mesh]
file = '/home/lindsayad/gdrive/gmsh-scripts/2d-msr.msh'
[]
[Variables]
[group1]
order = FIRST
family = LAGRANGE
[]
[group2]
order = FIRST
family = LAGRANGE
[]
[temp]
order = FIRST
family = LAGRANGE
scaling = 1e-6
[]
[]
[Kernels]
# Neutronics
[diff_group1]
type = GroupDiffusion
variable = group1
group_number = 1
[]
[diff_group2]
type = GroupDiffusion
variable = group2
group_number = 2
[]
[sigma_r_group1]
type = SigmaR
variable = group1
group_number = 1
[]
[sigma_r_group2]
type = SigmaR
variable = group2
group_number = 2
[]
[inscatter_group1]
type = InScatter
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[]
[inscatter_group2]
type = InScatter
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[]
[fission_source_group1]
type = CoupledFissionEigenKernel
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[]
[fission_source_group2]
type = CoupledFissionEigenKernel
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[]
# Temperature
[temp_cond]
type = MatDiffusion
variable = temp
prop_name = 'k'
save_in = 'diffus_resid tot_resid'
[]
# [./temp_cond]
# type = HeatConduction
# diffusion_coefficient_name = k
# diffusion_coefficient_dT_name = d_k_d_temp
# use_displaced_mesh = false
# variable = temp
# [../]
[temp_source]
type = FissionHeatSource
tot_fissions = tot_fissions
variable = temp
save_in = 'src_resid tot_resid'
[]
[]
[AuxVariables]
[Qf]
family = MONOMIAL
order = CONSTANT
[]
[diffus_temp]
family = MONOMIAL
order = CONSTANT
[]
[diffus_resid]
family = LAGRANGE
order = FIRST
[]
[src_resid]
family = LAGRANGE
order = FIRST
[]
[bc_resid]
family = LAGRANGE
order = FIRST
[]
[tot_resid]
family = LAGRANGE
order = FIRST
[]
[]
[AuxKernels]
[Qf]
type = FissionHeatSourceAux
variable = Qf
tot_fissions = tot_fissions
[]
[diffus_temp]
type = MatDiffusionAux
variable = diffus_temp
diffuse_var = temp
prop_name = 'k'
[]
[]
[Materials]
[fuel]
type = GenericMoltresMaterial
block = 'fuel'
property_tables_root = '/home/lindsayad/serpent/core/examples/serpent-input/msre/msr2g_enrU_mod_953_fuel_interp_'
num_groups = 2
prop_names = 'k'
prop_values = '.0123' # Cammi 2011 at 908 K
# prop_names = 'k d_k_d_temp'
# prop_values = '.0123 0' # Cammi 2011 at 908 K
[]
[moder]
type = GenericMoltresMaterial
block = 'moder'
property_tables_root = '/home/lindsayad/serpent/core/examples/serpent-input/msre/msr2g_enrU_fuel_922_mod_interp_'
num_groups = 2
prop_names = 'k'
prop_values = '.312' # Cammi 2011 at 908 K
# prop_names = 'k d_k_d_temp'
# prop_values = '.312 0' # Cammi 2011 at 908 K
[]
[]
[BCs]
[temp]
boundary = 'boundary'
type = DirichletBC
variable = temp
value = 900
save_in = 'bc_resid tot_resid'
[]
# [./temp]
# boundary = boundary
# type = VacuumBC
# variable = temp
# [../]
[]
[Executioner]
# type = NonlinearEigen
# free_power_iterations = 4
# source_abs_tol = 1e-12
# source_rel_tol = 1e-8
# output_after_power_iterations = true
type = InversePowerMethod
max_power_iterations = 50
xdiff = 'group1diff'
bx_norm = 'bnorm'
k0 = 1.0
pfactor = 1e-2
l_max_its = 100
solve_type = 'PJFNK'
# solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
# petsc_options_iname = '-pc_type -sub_pc_type'
# petsc_options_value = 'asm lu'
[]
[Preconditioning]
[SMP]
type = SMP
# full = true
[]
[]
[Postprocessors]
[bnorm]
type = ElmIntegTotFissNtsPostprocessor
group_fluxes = 'group1 group2'
execute_on = linear
[]
[tot_fissions]
type = ElmIntegTotFissPostprocessor
execute_on = linear
[]
[group1norm]
type = ElementIntegralVariablePostprocessor
variable = group1
execute_on = linear
[]
[group2norm]
type = ElementIntegralVariablePostprocessor
variable = group2
execute_on = linear
[]
[group1max]
type = NodalExtremeValue
value_type = max
variable = group1
execute_on = timestep_end
[]
[group2max]
type = NodalExtremeValue
value_type = max
variable = group2
execute_on = timestep_end
[]
[src_resid_post]
type = NodalL2Norm
variable = src_resid
execute_on = nonlinear
[]
[diffus_resid_post]
type = NodalL2Norm
variable = diffus_resid
execute_on = nonlinear
[]
[bc_resid_post]
type = NodalL2Norm
variable = bc_resid
execute_on = nonlinear
[]
[tot_resid_post]
type = NodalL2Norm
variable = tot_resid
execute_on = nonlinear
[]
[group1diff]
type = ElementL2Diff
variable = group1
execute_on = 'linear timestep_end'
use_displaced_mesh = false
[]
[]
[Outputs]
[out]
type = Exodus
execute_on = 'timestep_end'
[]
[]
[Debug]
show_var_residual_norms = true
[]
[ICs]
[temp_ic]
type = ConstantIC
variable = temp
value = 900
[]
[]
(problems/debug-precursor-action.i)
flow_velocity = 147 # Cammi 147 cm/s
[GlobalParams]
num_groups = 2
num_precursor_groups = 1
group_fluxes = 'group1 group2'
# MSRE full power = 10 MW; core volume 90 ft3
power = 200000
[]
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
block_id = 0
block_name = 'fuel'
[]
[MeshModifiers]
[mod]
type = SubdomainBoundingBox
bottom_left = '0.5 0 0'
top_right = '1 1 0'
block_id = 1
block_name = 'moder'
[]
[]
[Variables]
# [./pre1]
# [../]
[group1]
order = FIRST
family = LAGRANGE
[]
[group2]
order = FIRST
family = LAGRANGE
[]
[temp]
order = FIRST
family = LAGRANGE
scaling = 1e-3
[]
[]
[Kernels]
# Neutronics
# [./time_group1]
# type = NtTimeDerivative
# grou_number = 1
# variable = group1
# [../]
# [./time_group2]
# type = NtTimeDerivative
# grou_number = 2
# variable = group2
# [../]
[diff_group1]
type = GroupDiffusion
variable = group1
group_number = 1
temperature = temp
[]
[diff_group2]
type = GroupDiffusion
variable = group2
group_number = 2
temperature = temp
[]
[sigma_r_group1]
type = SigmaR
variable = group1
group_number = 1
temperature = temp
[]
[sigma_r_group2]
type = SigmaR
variable = group2
group_number = 2
temperature = temp
[]
[inscatter_group1]
type = InScatter
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
temperature = temp
[]
[inscatter_group2]
type = InScatter
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
temperature = temp
[]
[fission_source_group1]
type = CoupledFissionEigenKernel
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
temperature = temp
[]
[fission_source_group2]
type = CoupledFissionEigenKernel
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
temperature = temp
[]
# [./fission_source_group1]
# type = CoupledFissionKernel
# variable = group1
# group_number = 1
# num_groups = 2
# group_fluxes = 'group1 group2'
# temperature = temp
# [../]
# [./fission_source_group2]
# type = CoupledFissionKernel
# variable = group2
# group_number = 2
# num_groups = 2
# group_fluxes = 'group1 group2'
# temperature = temp
# [../]
# Temperature
# [./temp_cond]
# type = MatDiffusion
# variable = temp
# prop_name = 'k'
# [../]
[temp_flow_fuel]
block = 'fuel'
type = MatINSTemperatureRZ
variable = temp
rho = 'rho'
k = 'k'
cp = 'cp'
uz = ${flow_velocity}
[]
[temp_flow_moder]
block = 'moder'
type = MatINSTemperatureRZ
variable = temp
rho = 'rho'
k = 'k'
cp = 'cp'
[]
[temp_source]
type = FissionHeatSource
tot_fissions = tot_fissions
variable = temp
[]
# # Delayed neutron precursors
# [./pre1_source]
# type = PrecursorSource
# variable = pre1
# precursor_group_number = 1
# temperature = temp
# [../]
# [./pre1_decay]
# type = PrecursorDecay
# variable = pre1
# precursor_group_number = 2
# temperature = temp
# [../]
[]
[Precursors]
var_name_base = pre
# v_def = ${flow_velocity}
# block = 'fuel'
inlet_boundary = 'bottom'
inlet_boundary_condition = 'DirichletBC'
inlet_dirichlet_value = 0
outlet_boundary = 'top'
T = temp
[]
[AuxVariables]
[Qf]
family = MONOMIAL
order = CONSTANT
[]
# [./diffus_temp]
# family = MONOMIAL
# order = CONSTANT
# [../]
# [./diffus_resid]
# family = LAGRANGE
# order = FIRST
# [../]
# [./src_resid]
# family = LAGRANGE
# order = FIRST
# [../]
# [./bc_resid]
# family = LAGRANGE
# order = FIRST
# [../]
# [./tot_resid]
# family = LAGRANGE
# order = FIRST
# [../]
[]
[AuxKernels]
[Qf]
type = FissionHeatSourceAux
variable = Qf
tot_fissions = tot_fissions
[]
# [./diffus_temp]
# type = MatDiffusionAux
# variable = diffus_temp
# diffuse_var = temp
# prop_name = 'k'
# [../]
[]
[Materials]
[fuel]
type = GenericMoltresMaterial
block = 'fuel'
property_tables_root = 'msr2g_enrU_mod_953_fuel_interp_'
num_groups = 2
num_precursor_groups = 8
prop_names = 'k rho cp'
prop_values = '.0123 3.327e-3 1357' # Cammi 2011 at 908 K
temperature = temp
[]
[moder]
type = GenericMoltresMaterial
block = 'moder'
property_tables_root = 'msr2g_enrU_fuel_922_mod_interp_'
num_groups = 2
num_precursor_groups = 8
prop_names = 'k rho cp'
prop_values = '.312 1.843e-3 1760' # Cammi 2011 at 908 K
temperature = temp
[]
[]
[BCs]
[temp_inlet]
boundary = 'bottom'
type = DirichletBC
variable = temp
value = 900
[]
[temp_outlet]
boundary = 'top'
type = MatINSTemperatureNoBCBC
variable = temp
k = 'k'
[]
[group1_vacuum]
type = VacuumBC
variable = group1
boundary = 'top bottom'
[]
[group2_vacuum]
type = VacuumBC
variable = group2
boundary = 'top bottom'
[]
[]
[Problem]
type = FEProblem
coord_type = RZ
[]
[Executioner]
type = NonlinearEigen
free_power_iterations = 2
source_abs_tol = 1e-12
source_rel_tol = 1e-8
output_after_power_iterations = false
# type = InversePowerMethod
# max_power_iterations = 50
# xdiff = 'group1diff'
bx_norm = 'bnorm'
k0 = 1.7
pfactor = 1e-2
l_max_its = 100
line_search = none
# solve_type = 'PJFNK'
solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor -snes_test_display'
# This system will not converge with default preconditioning; need to use asm
petsc_options_iname = '-pc_type -sub_pc_type -sub_ksp_type -pc_asm_overlap -ksp_gmres_restart'
petsc_options_value = 'asm lu preonly 2 31'
# petsc_options_iname = '-snes_type'
# petsc_options_value = 'test'
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Postprocessors]
[bnorm]
type = ElmIntegTotFissNtsPostprocessor
group_fluxes = 'group1 group2'
execute_on = linear
[]
[tot_fissions]
type = ElmIntegTotFissPostprocessor
execute_on = linear
[]
[group1norm]
type = ElementIntegralVariablePostprocessor
variable = group1
execute_on = linear
[]
[group2norm]
type = ElementIntegralVariablePostprocessor
variable = group2
execute_on = linear
[]
[group1max]
type = NodalExtremeValue
value_type = max
variable = group1
execute_on = timestep_end
[]
[group2max]
type = NodalExtremeValue
value_type = max
variable = group2
execute_on = timestep_end
[]
[group1diff]
type = ElementL2Diff
variable = group1
execute_on = 'linear timestep_end'
use_displaced_mesh = false
[]
[]
[Outputs]
[out]
type = Exodus
execute_on = 'timestep_end'
[]
[]
[Debug]
show_var_residual_norms = true
[]
[ICs]
[temp_ic]
type = ConstantIC
variable = temp
value = 900
[]
[]
(problems/axi-nts-temp-delayed.i)
flow_velocity = 147 # Cammi 147 cm/s
[GlobalParams]
num_groups = 2
num_precursor_groups = 8
group_fluxes = 'group1 group2'
# MSRE full power = 10 MW; core volume 90 ft3
power = 200000
[]
[Mesh]
file = 'axisymm_cylinder.msh'
[]
[Variables]
# [./pre1]
# [../]
[group1]
order = FIRST
family = LAGRANGE
[]
[group2]
order = FIRST
family = LAGRANGE
[]
[temp]
order = FIRST
family = LAGRANGE
scaling = 1e-3
[]
[]
[Kernels]
# Neutronics
# [./time_group1]
# type = NtTimeDerivative
# grou_number = 1
# variable = group1
# [../]
# [./time_group2]
# type = NtTimeDerivative
# grou_number = 2
# variable = group2
# [../]
[diff_group1]
type = GroupDiffusion
variable = group1
group_number = 1
temperature = temp
[]
[diff_group2]
type = GroupDiffusion
variable = group2
group_number = 2
temperature = temp
[]
[sigma_r_group1]
type = SigmaR
variable = group1
group_number = 1
temperature = temp
[]
[sigma_r_group2]
type = SigmaR
variable = group2
group_number = 2
temperature = temp
[]
[inscatter_group1]
type = InScatter
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
temperature = temp
[]
[inscatter_group2]
type = InScatter
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
temperature = temp
[]
[fission_source_group1]
type = CoupledFissionEigenKernel
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
temperature = temp
[]
[fission_source_group2]
type = CoupledFissionEigenKernel
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
temperature = temp
[]
# [./fission_source_group1]
# type = CoupledFissionKernel
# variable = group1
# group_number = 1
# num_groups = 2
# group_fluxes = 'group1 group2'
# temperature = temp
# [../]
# [./fission_source_group2]
# type = CoupledFissionKernel
# variable = group2
# group_number = 2
# num_groups = 2
# group_fluxes = 'group1 group2'
# temperature = temp
# [../]
# Temperature
# [./temp_cond]
# type = MatDiffusion
# variable = temp
# prop_name = 'k'
# [../]
[temp_flow_fuel]
block = 'fuel'
type = MatINSTemperatureRZ
variable = temp
rho = 'rho'
k = 'k'
cp = 'cp'
uz = ${flow_velocity}
[]
[temp_flow_moder]
block = 'moder'
type = MatINSTemperatureRZ
variable = temp
rho = 'rho'
k = 'k'
cp = 'cp'
[]
[temp_source]
type = FissionHeatSource
tot_fissions = tot_fissions
variable = temp
[]
[]
[Precursors]
var_name_base = pre
v_def = ${flow_velocity}
block = 'fuel'
inlet_boundary = 'fuel_bottom'
inlet_boundary_condition = 'DirichletBC'
inlet_dirichlet_value = -20
outlet_boundary = 'fuel_top'
T = temp
incompressible_flow = false
use_exp_form = true
initial_condition = -20
[]
[AuxVariables]
[Qf]
family = MONOMIAL
order = CONSTANT
[]
[pre1_lin]
family = MONOMIAL
order = CONSTANT
[]
# [./diffus_temp]
# family = MONOMIAL
# order = CONSTANT
# [../]
# [./diffus_resid]
# family = LAGRANGE
# order = FIRST
# [../]
# [./src_resid]
# family = LAGRANGE
# order = FIRST
# [../]
# [./bc_resid]
# family = LAGRANGE
# order = FIRST
# [../]
# [./tot_resid]
# family = LAGRANGE
# order = FIRST
# [../]
[]
[AuxKernels]
[Qf]
type = FissionHeatSourceAux
variable = Qf
tot_fissions = tot_fissions
[]
[pre1_lin]
variable = pre1_lin
density_log = pre1
type = Density
[]
# [./diffus_temp]
# type = MatDiffusionAux
# variable = diffus_temp
# diffuse_var = temp
# prop_name = 'k'
# [../]
[]
[Materials]
[fuel]
type = GenericMoltresMaterial
block = 'fuel'
property_tables_root = 'msr2g_enrU_mod_953_fuel_interp_'
num_groups = 2
num_precursor_groups = 8
prop_names = 'k rho cp'
prop_values = '.0123 3.327e-3 1357' # Cammi 2011 at 908 K
temperature = temp
[]
[moder]
type = GenericMoltresMaterial
block = 'moder'
property_tables_root = 'msr2g_enrU_fuel_922_mod_interp_'
num_groups = 2
num_precursor_groups = 8
prop_names = 'k rho cp'
prop_values = '.312 1.843e-3 1760' # Cammi 2011 at 908 K
temperature = temp
[]
[]
[BCs]
[temp_inlet]
boundary = 'fuel_bottom graphite_bottom'
type = DirichletBC
variable = temp
value = 900
[]
[temp_outlet]
boundary = 'fuel_top'
type = MatINSTemperatureNoBCBC
variable = temp
k = 'k'
[]
[group1_vacuum]
type = VacuumBC
variable = group1
boundary = 'fuel_top graphite_top fuel_bottom graphite_bottom'
[]
[group2_vacuum]
type = VacuumBC
variable = group2
boundary = 'fuel_top graphite_top fuel_bottom graphite_bottom'
[]
[]
[Problem]
type = FEProblem
coord_type = RZ
[]
[Executioner]
type = NonlinearEigen
free_power_iterations = 2
source_abs_tol = 1e-12
source_rel_tol = 1e-8
output_after_power_iterations = false
# type = InversePowerMethod
# max_power_iterations = 50
# xdiff = 'group1diff'
bx_norm = 'bnorm'
k0 = 1.7
pfactor = 1e-2
l_max_its = 100
# line_search = none
solve_type = 'PJFNK'
# solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
# This system will not converge with default preconditioning; need to use asm
petsc_options_iname = '-pc_type -sub_pc_type -sub_ksp_type -pc_asm_overlap -ksp_gmres_restart'
petsc_options_value = 'asm lu preonly 2 31'
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Postprocessors]
[bnorm]
type = ElmIntegTotFissNtsPostprocessor
group_fluxes = 'group1 group2'
execute_on = linear
[]
[tot_fissions]
type = ElmIntegTotFissPostprocessor
execute_on = linear
[]
[group1norm]
type = ElementIntegralVariablePostprocessor
variable = group1
execute_on = linear
[]
[group2norm]
type = ElementIntegralVariablePostprocessor
variable = group2
execute_on = linear
[]
[group1max]
type = NodalExtremeValue
value_type = max
variable = group1
execute_on = timestep_end
[]
[group2max]
type = NodalExtremeValue
value_type = max
variable = group2
execute_on = timestep_end
[]
[group1diff]
type = ElementL2Diff
variable = group1
execute_on = 'linear timestep_end'
use_displaced_mesh = false
[]
[]
[Outputs]
[out]
type = Exodus
execute_on = 'timestep_end'
[]
[]
[Debug]
show_var_residual_norms = true
[]
[ICs]
[temp_ic]
type = ConstantIC
variable = temp
value = 900
[]
[]
(problems/msr-couple-temp-nts-2x2.i)
[GlobalParams]
num_groups = 2
temperature = temp
group_fluxes = 'group1 group2'
[../]
[Mesh]
file = '/home/lindsayad/gdrive/gmsh-scripts/msr-small.msh'
[../]
[Variables]
[./group1]
order = FIRST
family = LAGRANGE
[../]
[./group2]
order = FIRST
family = LAGRANGE
[../]
[./temp]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
# Neutronics
[./diff_group1]
type = GroupDiffusion
variable = group1
group_number = 1
[../]
[./diff_group2]
type = GroupDiffusion
variable = group2
group_number = 2
[../]
[./sigma_r_group1]
type = SigmaR
variable = group1
group_number = 1
[../]
[./sigma_r_group2]
type = SigmaR
variable = group2
group_number = 2
[../]
[./inscatter_group1]
type = InScatter
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[../]
[./inscatter_group2]
type = InScatter
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[../]
[./fission_source_group1]
type = CoupledFissionEigenKernel
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[../]
[./fission_source_group2]
type = CoupledFissionEigenKernel
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[../]
# Temperature
[./temp_cond]
type = HeatConduction
diffusion_coefficient_name = k
diffusion_coefficient_dT_name = d_k_d_temp
use_displaced_mesh = false
variable = temp
[../]
[./temp_source]
type = FissionHeatSource
tot_fissions = tot_fissions
# MSRE full power = 10 MW; core volume 90 ft3
power = 7000 # Watts
variable = temp
[../]
[]
[Materials]
[./fuel]
type = GenericMoltresMaterial
block = 'fuel'
property_tables_root = '/home/lindsayad/serpent/core/examples/serpent-input/msre/msr2g_enrU_mod_953_fuel_interp_'
num_groups = 2
prop_names = 'k d_k_d_temp'
prop_values = '.0123 0' # Cammi 2011 at 908 K
[../]
[./moder]
type = GenericMoltresMaterial
block = 'moder'
property_tables_root = '/home/lindsayad/serpent/core/examples/serpent-input/msre/msr2g_enrU_fuel_922_mod_interp_'
num_groups = 2
prop_names = 'k d_k_d_temp'
prop_values = '.312 0' # Cammi 2011 at 908 K
[../]
[]
[BCs]
[./temp]
boundary = boundary
type = DirichletBC
variable = temp
value = 900
[../]
# [./temp]
# boundary = boundary
# type = VacuumBC
# variable = temp
# [../]
[]
[Executioner]
type = NonlinearEigen
bx_norm = 'bnorm'
free_power_iterations = 2
source_abs_tol = 1e-12
source_rel_tol = 1e-50
k0 = 1.0
output_after_power_iterations = true
# solve_type = 'PJFNK'
solve_type = 'NEWTON'
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
[../]
[]
[Postprocessors]
[./bnorm]
type = ElmIntegTotFissNtsPostprocessor
group_fluxes = 'group1 group2'
execute_on = linear
[../]
[./tot_fissions]
type = ElmIntegTotFissPostprocessor
execute_on = linear
[../]
[]
[Outputs]
[./out]
type = Exodus
execute_on = 'timestep_end'
[../]
[]
(problems/gen-mesh-one-material.i)
[GlobalParams]
num_groups = 2
temperature = temp
group_fluxes = 'group1 group2'
# MSRE full power = 10 MW; core volume 90 ft3
power = 10
[]
[Mesh]
type = GeneratedMesh
dim = 2
xmax = 6
ymax = 6
nx = 15
ny = 15
[]
[Variables]
[group1]
order = FIRST
family = LAGRANGE
[]
[group2]
order = FIRST
family = LAGRANGE
[]
[temp]
order = FIRST
family = LAGRANGE
scaling = 1e-6
[]
[]
[Kernels]
# Neutronics
[diff_group1]
type = GroupDiffusion
variable = group1
group_number = 1
[]
[diff_group2]
type = GroupDiffusion
variable = group2
group_number = 2
[]
[sigma_r_group1]
type = SigmaR
variable = group1
group_number = 1
[]
[sigma_r_group2]
type = SigmaR
variable = group2
group_number = 2
[]
[inscatter_group1]
type = InScatter
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[]
[inscatter_group2]
type = InScatter
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[]
[fission_source_group1]
type = CoupledFissionEigenKernel
variable = group1
group_number = 1
num_groups = 2
group_fluxes = 'group1 group2'
[]
[fission_source_group2]
type = CoupledFissionEigenKernel
variable = group2
group_number = 2
num_groups = 2
group_fluxes = 'group1 group2'
[]
# Temperature
[temp_cond]
type = MatDiffusion
variable = temp
prop_name = 'k'
save_in = 'diffus_resid tot_resid'
[]
# [./temp_cond]
# type = HeatConduction
# diffusion_coefficient_name = k
# diffusion_coefficient_dT_name = d_k_d_temp
# use_displaced_mesh = false
# variable = temp
# [../]
[temp_source]
type = FissionHeatSource
tot_fissions = tot_fissions
variable = temp
save_in = 'src_resid tot_resid'
[]
[]
[AuxVariables]
[Qf]
family = MONOMIAL
order = CONSTANT
[]
[diffus_temp]
family = MONOMIAL
order = CONSTANT
[]
[diffus_resid]
family = LAGRANGE
order = FIRST
[]
[src_resid]
family = LAGRANGE
order = FIRST
[]
[bc_resid]
family = LAGRANGE
order = FIRST
[]
[tot_resid]
family = LAGRANGE
order = FIRST
[]
[]
[AuxKernels]
[Qf]
type = FissionHeatSourceAux
variable = Qf
tot_fissions = tot_fissions
[]
[diffus_temp]
type = MatDiffusionAux
variable = diffus_temp
diffuse_var = temp
prop_name = 'k'
[]
[]
[Materials]
[fuel]
type = GenericMoltresMaterial
block = 0
property_tables_root = '/home/lindsayad/serpent/core/examples/serpent-input/msre/msr2g_enrU_mod_953_fuel_interp_'
num_groups = 2
prop_names = 'k'
prop_values = '.0123' # Cammi 2011 at 908 K
# prop_names = 'k d_k_d_temp'
# prop_values = '.0123 0' # Cammi 2011 at 908 K
[]
[]
[BCs]
[temp]
boundary = 'left right top bottom'
type = DirichletBC
variable = temp
value = 900
save_in = 'bc_resid tot_resid'
[]
# [./temp]
# boundary = boundary
# type = VacuumBC
# variable = temp
# [../]
[]
[Executioner]
# type = NonlinearEigen
# free_power_iterations = 4
# source_abs_tol = 1e-12
# source_rel_tol = 1e-8
# output_after_power_iterations = true
type = InversePowerMethod
max_power_iterations = 50
xdiff = 'group1diff'
bx_norm = 'bnorm'
k0 = 1.0
pfactor = 1e-2
l_max_its = 100
solve_type = 'PJFNK'
# solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
# petsc_options_iname = '-pc_type -sub_pc_type'
# petsc_options_value = 'asm lu'
[]
[Preconditioning]
[SMP]
type = SMP
# full = true
[]
[]
[Postprocessors]
[bnorm]
type = ElmIntegTotFissNtsPostprocessor
group_fluxes = 'group1 group2'
execute_on = linear
[]
[tot_fissions]
type = ElmIntegTotFissPostprocessor
execute_on = linear
[]
[group1norm]
type = ElementIntegralVariablePostprocessor
variable = group1
execute_on = linear
[]
[group2norm]
type = ElementIntegralVariablePostprocessor
variable = group2
execute_on = linear
[]
[group1max]
type = NodalExtremeValue
value_type = max
variable = group1
execute_on = timestep_end
[]
[group2max]
type = NodalExtremeValue
value_type = max
variable = group2
execute_on = timestep_end
[]
[src_resid_post]
type = NodalL2Norm
variable = src_resid
execute_on = nonlinear
[]
[diffus_resid_post]
type = NodalL2Norm
variable = diffus_resid
execute_on = nonlinear
[]
[bc_resid_post]
type = NodalL2Norm
variable = bc_resid
execute_on = nonlinear
[]
[tot_resid_post]
type = NodalL2Norm
variable = tot_resid
execute_on = nonlinear
[]
[group1diff]
type = ElementL2Diff
variable = group1
execute_on = 'linear timestep_end'
use_displaced_mesh = false
[]
[]
[Outputs]
[out]
type = Exodus
execute_on = 'timestep_end'
[]
[]
[Debug]
show_var_residual_norms = true
[]
[ICs]
[temp_ic]
type = ConstantIC
variable = temp
value = 900
[]
[]