HTR-10 Griffin Neutronics Model
Contact: Javier Ortensi, [email protected]
Model link: HTR-10 Griffin Model
The HTR-10 is a small pebble-bed test reactor rated at a thermal power of MWt intended as a steppingstone for the development of PBR technology in China. The HTR-10 is located at the Institute of Nuclear and New Energy Technology (INET) and achieved initial criticality on 1 December 2000. The design of the HTR-10 reactor represents the design features of the modular High-Temperature Gas-cooled Reactor (HTGR) which is primarily characterized by its inherent safety features. The reactor geometry is depicted in Figure 1.

Figure 1: HTR-10 geometry with vessel (IRPEEP, 2006).
The following discussion of the reactor geometry and the benchmark based thereon are taken from (IRPEEP, 2006). This reference discusses the initial critical core configuration of the HTR-10. Overall, the relevant core for the neutronics, reflector and shielding regions (everything inside of the boronated carbon bricks and carbon bricks in Figure 1) are 6.1 meters tall and have a radius of 1.9 meters; the core containing the pebble bed is split into the upper part containing the cylindrical core and the upper cone; and a lower part containing the lower cone and the discharge tube. The upper part of the core during the initial critical contains a mix of 16,890 fuel and dummy graphite pebbles with a ratio of 57 : 43, while the lower part of the core only contains dummy pebbles. The benchmark report assumes a uniform packing fraction of 0.61 throughout the whole core.
In addition to the initial critical configuration, Ref.(IRPEEP, 2003) describes the full core configuration that features a significantly larger number of pebbles in the upper core region. This full core configuration is attained by loading more pebbles starting from the initial critical configuration until the core is capable of being operated at full power. The volume of the fuel pebbles in the upper core region is estimated to be m; at a packing fraction of and a pebble radius of cm, this corresponds to pebbles in the upper core region.
Griffin Model
Griffin is used to model the neutronics of the HTR. This section will present and describe the input files for the initial critical and full cores.
Initial Critical Core
The complete input file for the initial critical core is shown below.
# ==============================================================================
# Model description
# Application : Griffin
# ------------------------------------------------------------------------------
# Idaho Falls, INL, May 31st, 2022
# Author(s): Javier Ortensi
# If using or referring to this model, please cite as explained in
# https://mooseframework.inl.gov/virtual_test_bed/citing.html
# ==============================================================================
# - HTR Critical Configuration GRIFFIN neutronics input
# - MainApp
# ==============================================================================
# - The Model has been built based on [1-2].
# ------------------------------------------------------------------------------
# [1] Benchmark Analysis of the HTR-10 with the MAMMOTH Reactor Physics
# Application, J. Ortensi et al.
# [2] Evaluation of high temperature gas cooled reactor performance: benchmark
# analysis related to initial testing of the httr and htr-10.
# Technical Report IAEA-TECDOC-1382, International Reactor Physics
# Experiment Evaluation Project, 2003.
# ==============================================================================
# MODEL PARAMETERS
# ==============================================================================
absorption_blocks = '10 11 12 13 14 15 16 17 18 20 21 22 30 31 32 40 41 42 43 44 45 46 50 52'
fuel_blocks = '16 17 20 21 30 31'
non_fuel_blocks = '10 11 12 13 14 15 18 22 32 42 43 44 45 46 50 52'
TDC_blocks = '40 41'
# State ------------------------------------------------------------------------
# see state on xs library
# Power ------------------------------------------------------------------------
# ==============================================================================
# GLOBAL PARAMETERS
# ==============================================================================
[GlobalParams]
coupled_flux_groups = 'sflux_g0 sflux_g1 sflux_g2 sflux_g3 sflux_g4 sflux_g5 sflux_g6 sflux_g7 sflux_g8 sflux_g9'
scalar_flux = 'sflux_g0 sflux_g1 sflux_g2 sflux_g3 sflux_g4 sflux_g5 sflux_g6 sflux_g7 sflux_g8 sflux_g9'
library_file = '../data/xs/htr-10-XS.xml'
library_name = 'htr-10-critical'
grid_names = 'default'
grid = '1'
plus = 1
isotopes = 'pseudo'
densities = '1.0'
[]
# ==============================================================================
# GEOMETRY AND MESH
# ==============================================================================
[Mesh]
[fmg]
type = FileMeshGenerator
file = '../data/mesh/htr-10-critical-a-rev6.e'
exodus_extra_element_integers = 'eqv_id material_id'
[]
[eqvid]
type = ExtraElementIDCopyGenerator
input = fmg
source_extra_element_id = eqv_id
target_extra_element_ids = 'equivalence_id'
[]
# These modifiers are used to remove the boronated bricks that surround the core
# from the original mesh
[delete_bricks]
type = BlockDeletionGenerator
input = eqvid
block = '3 4 5'
[]
[sideset_side]
type = ParsedGenerateSideset
input = delete_bricks
combinatorial_geometry = 'sqrt(x*x+y*y) > 165.0'
new_sideset_name = 1
[]
[sideset_bot]
type = ParsedGenerateSideset
input = sideset_side
combinatorial_geometry = 'z < 41'
included_subdomains = '10 11 12'
normal = '0 0 -1'
new_sideset_name = 2
[]
[sideset_top]
type = ParsedGenerateSideset
input = sideset_bot
combinatorial_geometry = 'z > 490'
included_subdomains = '52 45 46'
normal = '0 0 1'
new_sideset_name = 3
[]
[]
[Equivalence]
type = SPH
transport_system = diff
compute_factors = false
equivalence_grid_names = 'default'
equivalence_grid = '1'
equivalence_library = '../data/sph/htr-10-Diff-SPH.xml'
[]
# ==============================================================================
# VARIABLES AND KERNELS
# ==============================================================================
[Variables]
[]
[Kernels]
[]
# ==============================================================================
# AUXVARIABLES AND AUXKERNELS
# ==============================================================================
[AuxVariables]
[AbsorptionRR]
order = FIRST
family = MONOMIAL
block = ${absorption_blocks}
[]
[NuFissionRR]
order = FIRST
family = MONOMIAL
block = ${fuel_blocks}
[]
[]
[AuxKernels]
[AbsorptionRR]
type = VectorReactionRate
block = ${absorption_blocks}
variable = AbsorptionRR
cross_section = sigma_absorption
scale_factor = power_scaling
dummies = _unscaled_total_power
execute_on = timestep_end
[]
[nuFissionRR]
type = VectorReactionRate
block = ${fuel_blocks}
variable = NuFissionRR
cross_section = nu_sigma_fission
scale_factor = power_scaling
dummies = _unscaled_total_power
execute_on = timestep_end
[]
[]
# ==============================================================================
# INITIAL CONDITIONS AND FUNCTIONS
# ==============================================================================
[ICs]
[]
[Functions]
[]
# ==============================================================================
# MATERIALS AND USER OBJECTS
# ==============================================================================
[Materials]
[fissile]
type = MixedMatIDNeutronicsMaterial
block = ${fuel_blocks}
[]
[non_fissile]
type = MixedMatIDNeutronicsMaterial
block = ${non_fuel_blocks}
[]
[TDC-materials]
type = MixedMatIDNeutronicsMaterial
block = ${TDC_blocks}
[]
[]
[PowerDensity]
power = 1
power_density_variable = power_density
integrated_power_postprocessor = integrated_power
power_scaling_postprocessor = power_scaling
family = MONOMIAL
order = CONSTANT
[]
# ==============================================================================
# BOUNDARY CONDITIONS
# ==============================================================================
[BCs]
[]
# ==============================================================================
# TRANSPORT SYSTEMS
# ==============================================================================
[TransportSystems]
particle = neutron
equation_type = eigenvalue
G = 10
VacuumBoundary = '1 2 3'
[diff]
scheme = CFEM-Diffusion
n_delay_groups = 0
family = LAGRANGE
order = FIRST
assemble_scattering_jacobian = true
assemble_fission_jacobian = true
[]
[]
# ==============================================================================
# EXECUTION PARAMETERS
# ==============================================================================
[Executioner]
type = Eigenvalue
# Preconditioned JFNK (default)
solve_type = 'PJFNKMO'
constant_matrices = true
petsc_options_iname = '-pc_type -pc_hypre_type -ksp_gmres_restart '
petsc_options_value = 'hypre boomeramg 100'
free_power_iterations = 1
nl_abs_tol = 1e-8
nl_rel_tol = 1e-7
[]
# ==============================================================================
# MULTIAPPS AND TRANSFER
# ==============================================================================
[MultiApps]
[]
[Transfers]
[]
# ==============================================================================
# RESTART
# ==============================================================================
[RestartVariables]
[]
# ==============================================================================
# POSTPROCESSORS DEBUG AND OUTPUTS
# ==============================================================================
[Postprocessors]
[power]
type = ElementIntegralVariablePostprocessor
variable = power_density
block = ${fuel_blocks}
[]
[nufission]
type = ElementIntegralVariablePostprocessor
variable = NuFissionRR
block = ${fuel_blocks}
[]
[absorption]
type = ElementIntegralVariablePostprocessor
variable = AbsorptionRR
block = ${absorption_blocks}
[]
[RR_Generation]
type = FluxRxnIntegral
cross_section = nu_sigma_fission
block = ${fuel_blocks}
[]
[RR_Absorption]
type = FluxRxnIntegral
cross_section = sigma_absorption
block = ${absorption_blocks}
[]
[RR_Leakage]
type = PartialSurfaceCurrent
boundary = '1 2 3'
transport_system = diff
[]
[]
[Outputs]
file_base = out_HTR-10_critical
time_step_interval = 1
exodus = true
csv = true
[console]
type = Console
[]
[]
(htgr/htr10/steady/htr-10-critical.i)In the following sections, we will discuss each of the input blocks.
Model Parameters
Model parameters are specified in this block. This could include the reactor state point (fuel temperature, coolant temperature/density, etc.) or reactor power level. The cross-section library contains only one state point, so it is not defined here.
Global Parameters
Global parameters are parameters that may be referenced throughout the input file. This block is user specific with the purpose to simplify repeated variable entries. However, be careful to not override a default parameter option in a later block without meaning to.
Regarding the multi-group cross section file, we define the library name and file location. In addition, the cross-section parameterization on the xml file is set. Macroscopic cross sections are always defined with the "pseudo" denomination and a number density value of 1.
If these definitions do not make sense now, we will clarify their uses and meanings in the subsequent blocks.
[GlobalParams]
coupled_flux_groups = 'sflux_g0 sflux_g1 sflux_g2 sflux_g3 sflux_g4 sflux_g5 sflux_g6 sflux_g7 sflux_g8 sflux_g9'
scalar_flux = 'sflux_g0 sflux_g1 sflux_g2 sflux_g3 sflux_g4 sflux_g5 sflux_g6 sflux_g7 sflux_g8 sflux_g9'
library_file = '../data/xs/htr-10-XS.xml'
library_name = 'htr-10-critical'
grid_names = 'default'
grid = '1'
plus = 1
isotopes = 'pseudo'
densities = '1.0'
[]
(htgr/htr10/steady/htr-10-critical.i)Geometry and Mesh
In this section, we will cover the mesh inputs. The full input block can be found below.
[Mesh]
[fmg]
type = FileMeshGenerator
file = '../data/mesh/htr-10-critical-a-rev6.e'
exodus_extra_element_integers = 'eqv_id material_id'
[]
[eqvid]
type = ExtraElementIDCopyGenerator
input = fmg
source_extra_element_id = eqv_id
target_extra_element_ids = 'equivalence_id'
[]
# These modifiers are used to remove the boronated bricks that surround the core
# from the original mesh
[delete_bricks]
type = BlockDeletionGenerator
input = eqvid
block = '3 4 5'
[]
[sideset_side]
type = ParsedGenerateSideset
input = delete_bricks
combinatorial_geometry = 'sqrt(x*x+y*y) > 165.0'
new_sideset_name = 1
[]
[sideset_bot]
type = ParsedGenerateSideset
input = sideset_side
combinatorial_geometry = 'z < 41'
included_subdomains = '10 11 12'
normal = '0 0 -1'
new_sideset_name = 2
[]
[sideset_top]
type = ParsedGenerateSideset
input = sideset_bot
combinatorial_geometry = 'z > 490'
included_subdomains = '52 45 46'
normal = '0 0 1'
new_sideset_name = 3
[]
[]
(htgr/htr10/steady/htr-10-critical.i)The HTR geometry is input into CUBIT, an external code developed at Sandia National Labs, via a CAD model to generate the computation mesh. An internal INL tool is used to post-process the mesh to ensure consistency with the Serpent generated multi-group cross section library (i.e., material regions and equivalent regions). The resulting output is an Exodus file that can be read by MOOSE using the FileMeshGenerator type. We need to directly specify the ids on the mesh using exodus_extra_element_integers which will then read in the material ids and equivalent ids (for the SPH factors). This is performed with the ExtraElementIDCopyGenerator type by defining the "eqv_id" as the source and "equivalence_id" as the target.
[Mesh]
[fmg]
type = FileMeshGenerator
file = '../data/mesh/htr-10-critical-a-rev6.e'
exodus_extra_element_integers = 'eqv_id material_id'
[]
[]
(htgr/htr10/steady/htr-10-critical.i)[Mesh]
[eqvid]
type = ExtraElementIDCopyGenerator
input = fmg
source_extra_element_id = eqv_id
target_extra_element_ids = 'equivalence_id'
[]
[]
(htgr/htr10/steady/htr-10-critical.i)Next, we need to make some modifications to the mesh. In particular, we want to delete non-neutronics relevant blocks on the original mesh (i.e., the boronated bricks). This can be performed in two steps: 1) delete the boronated bricks and 2) define the new boundary. The deletion of the bricks requires the use of the BlockDeletionGenerator type. Here, we select the blocks "3 4 5" that represent the boronated bricks in the original mesh. A simple way to identify the block ids is to open the mesh file HTR10-critical-c-rev2.e
in an Exodus supported visualization tool (i.e., ParaView). Next, we define the new top, bottom, and radial boundary. The radial boundary or [./sideset_side]
, in the model, is given the ParsedGenerateSideset type. This new boundary is defined by a function introduced in combinatorial_geometry and given a new name with new_sideset_name.
The top and bottom boundaries are defined in a similar manner. Again, the ParsedGenerateSideset type is used to define a new boundary with a function declared in combinatorial_geometry. The optional parameter, normal, specifies the normal vector on sides that are added to the new mesh. Also note that the input parameter is the sequential order of the named sub-blocks [./*]
.
[Mesh]
# These modifiers are used to remove the boronated bricks that surround the core
# from the original mesh
[delete_bricks]
type = BlockDeletionGenerator
input = eqvid
block = '3 4 5'
[]
[sideset_side]
type = ParsedGenerateSideset
input = delete_bricks
combinatorial_geometry = 'sqrt(x*x+y*y) > 165.0'
new_sideset_name = 1
[]
[sideset_bot]
type = ParsedGenerateSideset
input = sideset_side
combinatorial_geometry = 'z < 41'
included_subdomains = '10 11 12'
normal = '0 0 -1'
new_sideset_name = 2
[]
[sideset_top]
type = ParsedGenerateSideset
input = sideset_bot
combinatorial_geometry = 'z > 490'
included_subdomains = '52 45 46'
normal = '0 0 1'
new_sideset_name = 3
[]
[]
(htgr/htr10/steady/htr-10-critical.i)Equivalence
Equivalence factors for this model are hosted on LFS. Please refer to LFS instructions to download them.
The Equivalence theory block/action is used to compute or apply equivalence factors. In this case, Super homogenization factors are applied. Griffin also supports the use of discontinuity factors. The equivalence library is first defined with equivalance_library which is an xml file generated by the ISOXML pre-processing code. Both the equivalence_grid_names and the equivalence_grid are then specified with values from the equivalence library. These variables define the SPH factor parameterization. In this example, the cross-section library consists of only one state point, defined with the name "default" and grid "1". These variable names are on the equivalence library xml file. The "diff", short for diffusion, system will be defined later.
[Equivalence]
type = SPH
transport_system = diff
compute_factors = false
equivalence_grid_names = 'default'
equivalence_grid = '1'
equivalence_library = '../data/sph/htr-10-Diff-SPH.xml'
[]
(htgr/htr10/steady/htr-10-critical.i)AuxVariables and AuxKernels
AuxVariables are variables that can be derived from the solution variables (i.e., scalar flux). An AuxKernel is a procedure that uses the solution variable to compute the AuxVariable (i.e., reaction rate).
There are two AuxVariables that are defined in this model: the absorption reaction rate and the neutron production reaction rate.
[AuxVariables]
[AbsorptionRR]
order = FIRST
family = MONOMIAL
block = ${absorption_blocks}
[]
[NuFissionRR]
order = FIRST
family = MONOMIAL
block = ${fuel_blocks}
[]
[]
(htgr/htr10/steady/htr-10-critical.i)The AuxKernels are locally defined with the names [./AbsorptionRR]
and [./nuFissionRR]
, and are of the VectorReactionRate type. The AuxVariable that the kernel acts on is defined with variable. We also give it a cross_section and scale_factor. Lastly, we tell it to execute_on the end of a time step.
[AuxKernels]
[AbsorptionRR]
type = VectorReactionRate
block = ${absorption_blocks}
variable = AbsorptionRR
cross_section = sigma_absorption
scale_factor = power_scaling
dummies = _unscaled_total_power
execute_on = timestep_end
[]
[nuFissionRR]
type = VectorReactionRate
block = ${fuel_blocks}
variable = NuFissionRR
cross_section = nu_sigma_fission
scale_factor = power_scaling
dummies = _unscaled_total_power
execute_on = timestep_end
[]
[]
(htgr/htr10/steady/htr-10-critical.i)Materials
The cross sections for this model are hosted on LFS. Please refer to LFS instructions to download them.
Material cross sections are specified with the multi-group cross section library defined by library_file in the [GlobalParameters]
block. In this example, we define three types of materials: fissile, non-fissile, and TDC. TDC is short for "tensor diffusion coefficient". For each of these materials we specify the type as a MixedMatIDNeutronicsMaterial. This type uses neutronics properties based on mixed isotopes from the multi-group cross section library. In this example, we have a pseudo mixture and density that is defined in the [GlobalParameters]
block. The MixedMatIDNeutronicsMaterial type differs from MixedNeutronicsMaterial because it reads in the material ids directly from the mesh. We must specify which blocks are fissile materials with the block card. For each block id, we also need to identify the densities and isotopes. For convenience, these parameters were put in the [GlobalParameters]
block.
In a similar process, non-fissile materials are defined in the next sub-block. These could be any regions that have non-fissile materials such as a reflector, structural materials, etc. The last sub-block is for materials with a tensor diffusion coefficient (TDC).
[Materials]
[fissile]
type = MixedMatIDNeutronicsMaterial
block = ${fuel_blocks}
[]
[non_fissile]
type = MixedMatIDNeutronicsMaterial
block = ${non_fuel_blocks}
[]
[TDC-materials]
type = MixedMatIDNeutronicsMaterial
block = ${TDC_blocks}
[]
[]
(htgr/htr10/steady/htr-10-critical.i)The power density is not calculated implicitly and must be defined in the [PowerDensity]
block. This block provides the power density with a user specified power and the proper flux scaling factor in the postprocessor. The two required arguments are power and power_density_variable. The power is in units of Watts. The power_density_variable simply sets the variable name. Optional arguments include the post-processors to use. The integrated_power_postprocessor sets the desired name for the integrated power post-processor. Likewise, the power_scaling_postprocessor sets the desired name for the power scaling post-processor. Lastly, the family and order parameters are the polynomial representations of the power density corresponding to the underlying FEM.
[PowerDensity]
power = 1
power_density_variable = power_density
integrated_power_postprocessor = integrated_power
power_scaling_postprocessor = power_scaling
family = MONOMIAL
order = CONSTANT
[]
(htgr/htr10/steady/htr-10-critical.i)Transport System
The transport system is the heart of the Griffin simulation. We define the transport particle and eigenvalue equation with particle and equation_type, respectively. The number of energy groups is set to "10" with G. Vacuum boundary conditions are imposed on side sets "1 2 3". The [./diff]
sub-block defines the method to solve the problem. We use Griffin's diffusion method built on the MOOSE continuous finite element method solver. Finite element parameters are set with the family and order, which define the family of functions used to approximate the solution and polynomial order. The last two options, assemble_scattering_jacobian and assemble_fission_jacobian are required to use the "PJFNKMO" method defined later in the [Executioner]
. This tells the finite element solver to not update the material cross sections at each linear iteration.
[TransportSystems]
particle = neutron
equation_type = eigenvalue
G = 10
VacuumBoundary = '1 2 3'
[diff]
scheme = CFEM-Diffusion
n_delay_groups = 0
family = LAGRANGE
order = FIRST
assemble_scattering_jacobian = true
assemble_fission_jacobian = true
[]
[]
(htgr/htr10/steady/htr-10-critical.i)Executioner
The [Executioner]
block tells the solver what type of problem it needs to solve. Here, we select Eigenvalue as the executioner type which will solve the eigenvalue problem for the criticality and multi-group scalar flux (eigenvalue and eigenvectors). We also specify to solve with a Preconditioned Jacobian Free Newton Krylov Matrix Only method as shown by setting solve_type equal to "PJFNKMO". The Matrix Only method forces the solver to not update material properties (i.e., cross sections) in the linear iterations of the solve. To use this solver, the constant_matrices parameter must be set to "true". There are several optional arguments that may also be included. For example, we define a few PETSc options and non-linear solver tolerances.
[Executioner]
type = Eigenvalue
# Preconditioned JFNK (default)
solve_type = 'PJFNKMO'
constant_matrices = true
petsc_options_iname = '-pc_type -pc_hypre_type -ksp_gmres_restart '
petsc_options_value = 'hypre boomeramg 100'
free_power_iterations = 1
nl_abs_tol = 1e-8
nl_rel_tol = 1e-7
[]
(htgr/htr10/steady/htr-10-critical.i)Post-processors, Debug, and Outputs
The last blocks are for post-processors, debug options, and outputs. A post-processor can be thought of as a function to derive a variable of interest from the solution. For example, the power density. This is defined with the ElementIntegralVariablePostprocessor type and the power density variable we created earlier, power_density. The fissile blocks must be identified to compute the power density with block. Integral properties such as the rate of neutron production, rate of neutron absorption, and rate of neutron leakage are integrated over each FEM element in the next two sub-blocks. The remaining sub-blocks define the total generation with a FluxRxnIntegral type.
[Postprocessors]
[power]
type = ElementIntegralVariablePostprocessor
variable = power_density
block = ${fuel_blocks}
[]
[nufission]
type = ElementIntegralVariablePostprocessor
variable = NuFissionRR
block = ${fuel_blocks}
[]
[absorption]
type = ElementIntegralVariablePostprocessor
variable = AbsorptionRR
block = ${absorption_blocks}
[]
[RR_Generation]
type = FluxRxnIntegral
cross_section = nu_sigma_fission
block = ${fuel_blocks}
[]
[RR_Absorption]
type = FluxRxnIntegral
cross_section = sigma_absorption
block = ${absorption_blocks}
[]
[RR_Leakage]
type = PartialSurfaceCurrent
boundary = '1 2 3'
transport_system = diff
[]
[]
(htgr/htr10/steady/htr-10-critical.i)Finally, the output block sets the output files from the simulation. Two of the most common options include the exodus and csv file. The Exodus file can be viewed with the same software as the mesh, but now will show the solution and solution derived quantities such as the multi-group scalar flux and power distribution. The csv file stores a summary of the solution that includes the criticality. The perf_graph parameter is helpful to evaluate the computational run time.
[Outputs]
file_base = out_HTR-10_critical
time_step_interval = 1
exodus = true
csv = true
[console]
type = Console
[]
[]
(htgr/htr10/steady/htr-10-critical.i)Full Core
The complete input file for the full core is shown below.
# ==============================================================================
# Model description
# Application : Griffin
# ------------------------------------------------------------------------------
# Idaho Falls, INL, June 8th, 2022
# Author(s): Javier Ortensi
# If using or referring to this model, please cite as explained in
# https://mooseframework.inl.gov/virtual_test_bed/citing.html
# ==============================================================================
# - HTR Full Configuration GRIFFIN neutronics input
# - MainApp
# ==============================================================================
# - The Model has been built based on [1-2].
# ------------------------------------------------------------------------------
# [1] Benchmark Analysis of the HTR-10 with the MAMMOTH Reactor Physics
# Application, J. Ortensi et al.
# [2] Evaluation of high temperature gas cooled reactor performance: benchmark
# analysis related to initial testing of the httr and htr-10.
# Technical Report IAEA-TECDOC-1382, International Reactor Physics
# Experiment Evaluation Project, 2003.
# ==============================================================================
# MODEL PARAMETERS
# ==============================================================================
absorption_blocks = '10 11 12 13 14 15 16 17 18 20 21 22 30 31 32 40 41 42 43 44 45 46 50 52'
fuel_blocks = '16 17 20 21 30 31'
non_fuel_blocks = '10 11 12 13 14 15 18 22 32 40 41 42 43 44 45 46 50 52'
# State ------------------------------------------------------------------------
# see state on xs library
# Power ------------------------------------------------------------------------
# ==============================================================================
# GLOBAL PARAMETERS
# ==============================================================================
[GlobalParams]
coupled_flux_groups = 'sflux_g0 sflux_g1 sflux_g2 sflux_g3 sflux_g4 sflux_g5 sflux_g6 sflux_g7 sflux_g8 sflux_g9'
scalar_flux = 'sflux_g0 sflux_g1 sflux_g2 sflux_g3 sflux_g4 sflux_g5 sflux_g6 sflux_g7 sflux_g8 sflux_g9'
#Library name options include 'htr-10-full-1RI', 'htr-10-full-393K', 'htr-10-full-523K',
# 'htr-10-full-ARI', and 'htr-10-full-ARO'.
library_file = '../data/xs/htr-10-XS.xml'
library_name = 'htr-10-full-ARO'
grid = '1'
grid_names = 'default'
isotopes = 'pseudo'
densities = '1.0'
plus = 1
[]
# ==============================================================================
# GEOMETRY AND MESH
# ==============================================================================
[Mesh]
[fmg]
type = FileMeshGenerator
file = '../data/mesh/htr-10-full-a-rev3.e'
exodus_extra_element_integers = 'eqv_id material_id'
[]
[eqvid]
type = ExtraElementIDCopyGenerator
input = fmg
source_extra_element_id = eqv_id
target_extra_element_ids = 'equivalence_id'
[]
uniform_refine = 0
# These modifiers are used to remove the boronated bricks that surround the core
# from the original mesh
[delete_bricks]
type = BlockDeletionGenerator
input = eqvid
block = '3 4 5'
[]
[sideset_side]
type = ParsedGenerateSideset
input = delete_bricks
combinatorial_geometry = 'sqrt(x*x+y*y) > 165.0'
new_sideset_name = 1
[]
[sideset_bot]
type = ParsedGenerateSideset
input = sideset_side
combinatorial_geometry = 'z < 41'
included_subdomains = '10 11 12'
normal = '0 0 -1'
new_sideset_name = 2
[]
[sideset_top]
type = ParsedGenerateSideset
input = sideset_bot
combinatorial_geometry = 'z > 490'
included_subdomains = '52 45 46'
normal = '0 0 1'
new_sideset_name = 3
[]
[]
[Equivalence]
type = SPH
transport_system = ts
compute_factors = false
equivalence_grid_names = 'default'
equivalence_grid = '1'
equivalence_library = '../data/sph/htr-10-Diff-SPH.xml'
[]
# ==============================================================================
# VARIABLES AND KERNELS
# ==============================================================================
[Variables]
[]
[Kernels]
[]
# ==============================================================================
# AUXVARIABLES AND AUXKERNELS
# ==============================================================================
[AuxVariables]
[AbsorptionRR]
order = FIRST
family = MONOMIAL
block = ${absorption_blocks}
[]
[NuFissionRR]
order = FIRST
family = MONOMIAL
block = ${fuel_blocks}
[]
[]
[AuxKernels]
[AbsorptionRR]
type = VectorReactionRate
block = ${absorption_blocks}
variable = AbsorptionRR
cross_section = sigma_absorption
scale_factor = power_scaling
dummies = _unscaled_total_power
execute_on = timestep_end
[]
[nuFissionRR]
type = VectorReactionRate
block = ${fuel_blocks}
variable = NuFissionRR
cross_section = nu_sigma_fission
scale_factor = power_scaling
dummies = _unscaled_total_power
execute_on = timestep_end
[]
[]
# ==============================================================================
# INITIAL CONDITIONS AND FUNCTIONS
# ==============================================================================
[ICs]
[]
[Functions]
[]
# ==============================================================================
# MATERIALS AND USER OBJECTS
# ==============================================================================
[Materials]
[fissile]
type = MixedMatIDNeutronicsMaterial
block = ${fuel_blocks}
[]
[non_fissile]
type = MixedMatIDNeutronicsMaterial
block = ${non_fuel_blocks}
[]
[]
[PowerDensity]
power = 1
power_density_variable = power_density
integrated_power_postprocessor = integrated_power
power_scaling_postprocessor = power_scaling
family = MONOMIAL
order = CONSTANT
[]
# ==============================================================================
# BOUNDARY CONDITIONS
# ==============================================================================
[BCs]
[]
# ==============================================================================
# TRANSPORT SYSTEMS
# ==============================================================================
[TransportSystems]
particle = neutron
equation_type = eigenvalue
G = 10
VacuumBoundary = '1 2 3'
[ts]
scheme = CFEM-Diffusion
n_delay_groups = 0
family = LAGRANGE
order = FIRST
assemble_scattering_jacobian = true
assemble_fission_jacobian = true
[]
[]
# ==============================================================================
# EXECUTION PARAMETERS
# ==============================================================================
[Executioner]
type = Eigenvalue
# Preconditioned JFNK (default)
solve_type = 'PJFNKMO'
constant_matrices = true
petsc_options_iname = '-pc_type -pc_hypre_type -ksp_gmres_restart '
petsc_options_value = 'hypre boomeramg 100'
free_power_iterations = 1
[]
# ==============================================================================
# MULTIAPPS AND TRANSFER
# ==============================================================================
[MultiApps]
[]
[Transfers]
[]
# ==============================================================================
# RESTART
# ==============================================================================
[RestartVariables]
[]
# ==============================================================================
# POSTPROCESSORS DEBUG AND OUTPUTS
# ==============================================================================
[Postprocessors]
[power]
type = ElementIntegralVariablePostprocessor
variable = power_density
block = ${fuel_blocks}
[]
[nufission]
type = ElementIntegralVariablePostprocessor
variable = NuFissionRR
block = ${fuel_blocks}
[]
[absorption]
type = ElementIntegralVariablePostprocessor
variable = AbsorptionRR
block = ${absorption_blocks}
[]
[RR_Generation]
type = FluxRxnIntegral
cross_section = nu_sigma_fission
block = ${fuel_blocks}
[]
[RR_Absorption]
type = FluxRxnIntegral
cross_section = sigma_absorption
block = ${absorption_blocks}
[]
[RR_Leakage]
type = PartialSurfaceCurrent
boundary = '1 2 3'
transport_system = ts
[]
[]
[Outputs]
file_base = out-HTR-10-full-Diff-SPH-ARO
time_step_interval = 1
exodus = true
csv = true
[console]
type = Console
[]
[]
(htgr/htr10/steady/htr-10-full.i)The input file to the full core is similar to that of the initial critical core. Still, we will examine any features that were not discussed previously.
Global Parameters
In this model, there are several multi-group cross section libraries to choose from. Each library contains a unique state point at which Serpent generated homogenized cross sections. The available options are all rods in (ARI), all rods out (ARO), one rod in (1RI), and temperatures at 393K and 523K. These cross-section libraries are in the "xs" folder.
If equivalence factors are to be used, make sure to select a library in the "sph" folder that is consistent with the cross-section library.
[GlobalParams]
coupled_flux_groups = 'sflux_g0 sflux_g1 sflux_g2 sflux_g3 sflux_g4 sflux_g5 sflux_g6 sflux_g7 sflux_g8 sflux_g9'
scalar_flux = 'sflux_g0 sflux_g1 sflux_g2 sflux_g3 sflux_g4 sflux_g5 sflux_g6 sflux_g7 sflux_g8 sflux_g9'
#Library name options include 'htr-10-full-1RI', 'htr-10-full-393K', 'htr-10-full-523K',
# 'htr-10-full-ARI', and 'htr-10-full-ARO'.
library_file = '../data/xs/htr-10-XS.xml'
library_name = 'htr-10-full-ARO'
grid = '1'
grid_names = 'default'
isotopes = 'pseudo'
densities = '1.0'
plus = 1
[]
(htgr/htr10/steady/htr-10-full.i)Acknowledgments
This model description section is summarized from (Ortensi et al., 2018) and prepared for the VTB by Thomas Folk.
References
- IRPEEP.
Evaluation of high temperature gas cooled reactor performance:benchmark analysis related to initial testing of the httr and htr-10.
Technical Report IAEA-TECDOC-1382, International Reactor Physics Experiment Evaluation Project, 2003.[BibTeX]
- IRPEEP.
Evaluation of the initial critical configuration of the HTR-10 pebble-bed reactor.
Technical Report HTR10-GCR-RESR-001, International Reactor Physics Experiment Evaluation Project, 2006.[BibTeX]
- J. Ortensi, S. Schunert, Y. Wang, V. Laboure, F. Gleicher, and R. Martineau.
Benchmark Analysis of the HTR-10 with the MAMMOTH Reactor Physics Application.
Technical Report INL/EXT-18-45453, Idaho National Laboratory, 2018.[BibTeX]