Divertor Monoblock During Pulsed Operation
TMAP8 is used to model tritium transport in a divertor monoblock to elucidate the effects of pulsed operation (up to fifty 1600-second plasma discharge and cool-down cycles) on the tritium in-vessel inventory source term and ex-vessel release term (i.e., tritium retention and permeation) for safety analysis. This example reproduces the results presented in Shimada et al. (2024).
General description of the simulation case and corresponding input file
Introduction
In a magnetic confinement fusion system, the divertor components will be subjected to intense particle and heat fluxes from plasma, as well as 14-MeV neutrons stemming from deuterium-tritium (D-T) reactions, thus creating extremely high temperature and hydrogen concentration gradients across tens of millimeters. Furthermore, the thermomechanical properties of the divertor materials will evolve over time in response to displacement damage and gas/solid transmutations. Thus, predictive modeling of tritium transport in the divertor poses enormous challenges, as it requires simulation capabilities capable of (1) multi-material configurations, (2) high thermal cycles, (3) complex 2D and 3D geometries, (4) and microstructural evolution resulting from displacement damage and gas/solid transmutations. For example, the ITER divertor, which consists of tungsten (W) monoblocks bonded to cooling tubes made of copper-chromium-zirconium (CuCrZr) alloy, is designed to withstand high heat fluxes (∼10 MW ) during steady-state plasma operation from intense D-T plasma flux (∼1024 ). In this example, we simulate tritium retention and permeation as well as thermal transport in an ITER-like 2D W-Cu-CuCrZr monoblock and demonstrate the first three (out of four) aforementioned required simulation capabilities for tritium transport in the divertor using TMAP8, as presented in Shimada et al. (2024). This simulation was also based on the cases published in Hodille et al. (2021).
The sections below describe the simulation details and explain how they translate into the example TMAP8 input file. Not every part of the input file is explained here. If the interested reader has more questions about this example case and how to modify this input file to adapt it to a different case, feel free to reach out to the TMAP8 development team on the TMAP8 GitHub discussion page.
Divertor monoblock geometry and mesh generation
Figure 1 shows the geometry, mesh, and temperature distribution of the 2D monoblock employed in Shimada et al. (2024) and in the present example. The 2D monoblock consisted of three materials:
W monoblock ( (mm) > 8.5),
Cu interlayer (7.5 < (mm) < 8.5),
CuCrZr tube (6.0 < (mm) < 7.5), and
HO cooling ( (mm) < 6.0),
where in a Cartesian coordinate system, and the center point (, ) = (0, 0) is defined at the center of the CuCrZr tube.

Figure 1: 2D monoblock: (left) geometry and mesh; (right) temperature distribution. This corresponds to Fig. 1 in Ref. Shimada et al. (2024).
MOOSE is equipped with a set of user-friendly, built-in mesh generators for creating meshes based on simple geometries (e.g., a monoblock). The built-in mesh generator ConcentricCircleMeshGenerator is used to create meshes for the 2D monoblock geometry. Only the left-hand (−14.0 < (mm) < 0 & −14.0 < (mm) < 14.0) portion of the monoblock is used, such that the mesh size is halved by assuming symmetry at the plane. This example uses a total of 9,418 nodes, creating 28,402 degrees of freedom in the non-linear system. The finer mesh size and time step are required to compute the trapping-limited diffusion regime, as described in TMAP8 verification case ver-1d. However, a set of relatively low detrapping energies (< 0.85 eV) is used here to reduce the computing time via the relatively coarse mesh size.
In the input file, the mesh is defined as:
[Mesh<<<{"href": "../../syntax/Mesh/index.html"}>>>]
[ccmg]
type = ConcentricCircleMeshGenerator<<<{"description": "This ConcentricCircleMeshGenerator source code is to generate concentric circle meshes.", "href": "../../source/meshgenerators/ConcentricCircleMeshGenerator.html"}>>>
num_sectors<<<{"description": "num_sectors % 2 = 0, num_sectors > 0Number of azimuthal sectors in each quadrant'num_sectors' must be an even number."}>>> = 36
rings<<<{"description": "Number of rings in each circle or in the enclosing square"}>>> = '1 30 20 110'
radii<<<{"description": "Radii of major concentric circles"}>>> = '${units 6 mm -> m} ${units 7.5 mm -> m} ${units 8.5 mm -> m}'
has_outer_square<<<{"description": "It determines if meshes for a outer square are added to concentric circle meshes."}>>> = on
pitch<<<{"description": "The enclosing square can be added to the completed concentric circle mesh.Elements are quad meshes."}>>> = '${units 28 mm -> m}'
portion<<<{"description": "Control of which part of mesh is created"}>>> = left_half
preserve_volumes<<<{"description": "Volume of concentric circles can be preserved using this function."}>>> = false
smoothing_max_it<<<{"description": "Number of Laplacian smoothing iterations"}>>> = 3
[]
[ssbsg1]
type = SideSetsBetweenSubdomainsGenerator<<<{"description": "MeshGenerator that creates a sideset composed of the nodes located between two or more subdomains.", "href": "../../source/meshgenerators/SideSetsBetweenSubdomainsGenerator.html"}>>>
input<<<{"description": "The mesh we want to modify"}>>> = ccmg
primary_block<<<{"description": "The primary set of blocks for which to draw a sideset between"}>>> = '4' # W
paired_block<<<{"description": "The paired set of blocks for which to draw a sideset between"}>>> = '3' # Cu
new_boundary<<<{"description": "The list of boundary names to create on the supplied subdomain"}>>> = '4to3'
[]
[ssbsg2]
type = SideSetsBetweenSubdomainsGenerator<<<{"description": "MeshGenerator that creates a sideset composed of the nodes located between two or more subdomains.", "href": "../../source/meshgenerators/SideSetsBetweenSubdomainsGenerator.html"}>>>
input<<<{"description": "The mesh we want to modify"}>>> = ssbsg1
primary_block<<<{"description": "The primary set of blocks for which to draw a sideset between"}>>> = '3' # Cu
paired_block<<<{"description": "The paired set of blocks for which to draw a sideset between"}>>> = '4' # W
new_boundary<<<{"description": "The list of boundary names to create on the supplied subdomain"}>>> = '3to4'
[]
[ssbsg3]
type = SideSetsBetweenSubdomainsGenerator<<<{"description": "MeshGenerator that creates a sideset composed of the nodes located between two or more subdomains.", "href": "../../source/meshgenerators/SideSetsBetweenSubdomainsGenerator.html"}>>>
input<<<{"description": "The mesh we want to modify"}>>> = ssbsg2
primary_block<<<{"description": "The primary set of blocks for which to draw a sideset between"}>>> = '3' # Cu
paired_block<<<{"description": "The paired set of blocks for which to draw a sideset between"}>>> = '2' # CuCrZr
new_boundary<<<{"description": "The list of boundary names to create on the supplied subdomain"}>>> = '3to2'
[]
[ssbsg4]
type = SideSetsBetweenSubdomainsGenerator<<<{"description": "MeshGenerator that creates a sideset composed of the nodes located between two or more subdomains.", "href": "../../source/meshgenerators/SideSetsBetweenSubdomainsGenerator.html"}>>>
input<<<{"description": "The mesh we want to modify"}>>> = ssbsg3
primary_block<<<{"description": "The primary set of blocks for which to draw a sideset between"}>>> = '2' # CuCrZr
paired_block<<<{"description": "The paired set of blocks for which to draw a sideset between"}>>> = '3' # Cu
new_boundary<<<{"description": "The list of boundary names to create on the supplied subdomain"}>>> = '2to3'
[]
[ssbsg5]
type = SideSetsBetweenSubdomainsGenerator<<<{"description": "MeshGenerator that creates a sideset composed of the nodes located between two or more subdomains.", "href": "../../source/meshgenerators/SideSetsBetweenSubdomainsGenerator.html"}>>>
input<<<{"description": "The mesh we want to modify"}>>> = ssbsg4
primary_block<<<{"description": "The primary set of blocks for which to draw a sideset between"}>>> = '2' # CuCrZr
paired_block<<<{"description": "The paired set of blocks for which to draw a sideset between"}>>> = '1' # H2O
new_boundary<<<{"description": "The list of boundary names to create on the supplied subdomain"}>>> = '2to1'
[]
[bdg]
type = BlockDeletionGenerator<<<{"description": "Mesh generator which removes elements from the specified subdomains", "href": "../../source/meshgenerators/BlockDeletionGenerator.html"}>>>
input<<<{"description": "The mesh we want to modify"}>>> = ssbsg5
block<<<{"description": "The list of blocks to be deleted"}>>> = '1' # H2O
[]
[]
Nomenclature of variables and physical parameters
Table 1 lists the variables and physical parameters used in this example with their units.
Table 1: Nomenclature of variables and physical parameters used in this example.
Symbol | Variable or Physical Property | Unit |
---|---|---|
Concentration of solute (mobile) species | m | |
Concentration of trapped species | m | |
Temperature | K | |
Diffusivity of solute (mobile) species | m s | |
Solubility of solute (mobile) species | m Pa | |
Total concentration of species | m | |
Concentration of empty trapping sites | m | |
Concentration of trapping sites | m | |
Trapping rate coefficient, | s | |
Release rate coefficient | s | |
Pre-exponential factor, | s | |
Detrapping energy | eV | |
Atomic number density | m | |
Density | g m | |
Specific heat | J kg K | |
Thermal conductivity | W m K |
Variables and governing equations
To simulate tritium and thermal transport, we define two sets of partial differential equations (PDEs). First, the strong form of the mass conservation equation for solute (mobile) T atoms, , is written as:
(1)
We use three sets of mass conservation equations to calculate the behaviors of solute T atoms in three different materials (i.e., W, Cu, CuCrZr). Second, the strong form of the conservation of energy equation is written as:
(2)
Table 1 lists all the symbols and definitions used in Eq. (1) and Eq. (2).
The next step is to convert these two strong-form PDEs into their weak forms by multiplying with a test function, , and integrating over a domain, with surface and outward-facing normal vector . Using the divergence theorem, one can obtain the weak form of the mass conservation equation (Eq. (1)) as follows:
(3)
Similarly, the weak form of the conservation of energy equation can be written by multiplying Eq. (2) with a test function, , and integrating over a domain, with surface and outward-facing normal vector . Using the divergence theorem, one can obtain the weak form of the conservation of energy equation as follows:
(4)
Then, to solve for the PDEs and physical phenomena, we can select appropriate kernels and boundary conditions (BCs) from MOOSE’s extensive library. The following three subsections describe each kernel, BC, and numerical method utilized in the present work.
In the input file, the variables are defined as:
[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
[temperature]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = '${units 300 K}'
[]
######################### Variables for W (block = 4)
[C_mobile_W]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = '${units 1.0e-20 m^-3}'
block = 4
[]
[C_trapped_W]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = '${units 1.0e-15 m^-3}'
block = 4
[]
######################### Variables for Cu (block = 3)
[C_mobile_Cu]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = '${units 5.0e-17 m^-3}'
block = 3
[]
[C_trapped_Cu]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = '${units 1.0e-15 m^-3}'
block = 3
[]
######################### Variables for CuCrZr (block = 2)
[C_mobile_CuCrZr]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = '${units 1.0e-15 m^-3}'
block = 2
[]
[C_trapped_CuCrZr]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = '${units 1.0e-15 m^-3}'
block = 2
[]
[]
Note the usage of initial_condition
and block
parameters in order to set initial conditions and material block applicability for each variable.
Kernels
For the mass conservation equation, we use ADTimeDerivative and ADMatDiffusion to represent the 1 and 3 terms of Eq. (3). The TMAP8 kernels TrappingNodalKernel and ReleasingNodalKernel represent the 4 and 5 terms of Eq. (3), respectively, and simulate the trapping/release behavior of hydrogen isotopes in/from trap sites. For the conservation of energy equation, SpecificHeatConductionTimeDerivative and HeatConduction solve the 1 and 3 terms of Eq. (4). ADTimeDerivative
and ADMatDiffusion
are kernels from the core MOOSE framework, and SpecificHeatConductionTimeDerivative
and HeatConduction
are kernels from the MOOSE Heat Transfer Module. These pre-made kernels are commonly re-used to represent time derivative, diffusion, and heat conduction terms in a given material within a variety of MOOSE-based simulations.
[Kernels<<<{"href": "../../syntax/Kernels/index.html"}>>>]
############################## Kernels for W (block = 4)
[diff_W]
type = ADMatDiffusion<<<{"description": "Diffusion equation kernel that takes an isotropic diffusivity from a material property", "href": "../../source/kernels/ADMatDiffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_W
diffusivity<<<{"description": "The diffusivity value or material property"}>>> = diffusivity_W
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[time_diff_W]
type = ADTimeDerivative<<<{"description": "The time derivative operator with the weak form of $(\\psi_i, \\frac{\\partial u_h}{\\partial t})$.", "href": "../../source/kernels/ADTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_W
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[coupled_time_W]
type = ScaledCoupledTimeDerivative<<<{"description": "Time derivative Kernel that acts on a coupled variable. Weak form: $(\\psi_i, \\frac{\\partial v_h}{\\partial t})$.", "href": "../../source/kernels/ScaledCoupledTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_W
v<<<{"description": "Coupled variable"}>>> = C_trapped_W
factor<<<{"description": "The factor by which to scale"}>>> = 1e0
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[heat_conduction_W]
type = HeatConduction<<<{"description": "Diffusive heat conduction term $-\\nabla\\cdot(k\\nabla T)$ of the thermal energy conservation equation", "href": "../../source/kernels/HeatConduction.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
diffusion_coefficient<<<{"description": "Property name of the diffusion coefficient"}>>> = thermal_conductivity_W
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[time_heat_conduction_W]
type = SpecificHeatConductionTimeDerivative<<<{"description": "Time derivative term $\\rho c_p \\frac{\\partial T}{\\partial t}$ of the heat equation with the specific heat $c_p$ and the density $\\rho$ as arguments.", "href": "../../source/kernels/SpecificHeatConductionTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
specific_heat<<<{"description": "Property name of the specific heat material property"}>>> = specific_heat_W
density<<<{"description": "Property name of the density material property"}>>> = density_W
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
############################## Kernels for Cu (block = 3)
[diff_Cu]
type = ADMatDiffusion<<<{"description": "Diffusion equation kernel that takes an isotropic diffusivity from a material property", "href": "../../source/kernels/ADMatDiffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_Cu
diffusivity<<<{"description": "The diffusivity value or material property"}>>> = diffusivity_Cu
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[time_diff_Cu]
type = ADTimeDerivative<<<{"description": "The time derivative operator with the weak form of $(\\psi_i, \\frac{\\partial u_h}{\\partial t})$.", "href": "../../source/kernels/ADTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_Cu
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[coupled_time_Cu]
type = ScaledCoupledTimeDerivative<<<{"description": "Time derivative Kernel that acts on a coupled variable. Weak form: $(\\psi_i, \\frac{\\partial v_h}{\\partial t})$.", "href": "../../source/kernels/ScaledCoupledTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_Cu
v<<<{"description": "Coupled variable"}>>> = C_trapped_Cu
factor<<<{"description": "The factor by which to scale"}>>> = 1e0
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[heat_conduction_Cu]
type = HeatConduction<<<{"description": "Diffusive heat conduction term $-\\nabla\\cdot(k\\nabla T)$ of the thermal energy conservation equation", "href": "../../source/kernels/HeatConduction.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
diffusion_coefficient<<<{"description": "Property name of the diffusion coefficient"}>>> = thermal_conductivity_Cu
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[time_heat_conduction_Cu]
type = SpecificHeatConductionTimeDerivative<<<{"description": "Time derivative term $\\rho c_p \\frac{\\partial T}{\\partial t}$ of the heat equation with the specific heat $c_p$ and the density $\\rho$ as arguments.", "href": "../../source/kernels/SpecificHeatConductionTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
specific_heat<<<{"description": "Property name of the specific heat material property"}>>> = specific_heat_Cu
density<<<{"description": "Property name of the density material property"}>>> = density_Cu
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
############################## Kernels for CuCrZr (block = 2)
[diff_CuCrZr]
type = ADMatDiffusion<<<{"description": "Diffusion equation kernel that takes an isotropic diffusivity from a material property", "href": "../../source/kernels/ADMatDiffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_CuCrZr
diffusivity<<<{"description": "The diffusivity value or material property"}>>> = diffusivity_CuCrZr
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[time_diff_CuCrZr]
type = ADTimeDerivative<<<{"description": "The time derivative operator with the weak form of $(\\psi_i, \\frac{\\partial u_h}{\\partial t})$.", "href": "../../source/kernels/ADTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_CuCrZr
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[coupled_time_CuCrZr]
type = ScaledCoupledTimeDerivative<<<{"description": "Time derivative Kernel that acts on a coupled variable. Weak form: $(\\psi_i, \\frac{\\partial v_h}{\\partial t})$.", "href": "../../source/kernels/ScaledCoupledTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_CuCrZr
v<<<{"description": "Coupled variable"}>>> = C_trapped_CuCrZr
factor<<<{"description": "The factor by which to scale"}>>> = 1e0
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[heat_conduction_CuCrZr]
type = HeatConduction<<<{"description": "Diffusive heat conduction term $-\\nabla\\cdot(k\\nabla T)$ of the thermal energy conservation equation", "href": "../../source/kernels/HeatConduction.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
diffusion_coefficient<<<{"description": "Property name of the diffusion coefficient"}>>> = thermal_conductivity_CuCrZr
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[time_heat_conduction_CuCrZr]
type = SpecificHeatConductionTimeDerivative<<<{"description": "Time derivative term $\\rho c_p \\frac{\\partial T}{\\partial t}$ of the heat equation with the specific heat $c_p$ and the density $\\rho$ as arguments.", "href": "../../source/kernels/SpecificHeatConductionTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
specific_heat<<<{"description": "Property name of the specific heat material property"}>>> = specific_heat_CuCrZr
density<<<{"description": "Property name of the density material property"}>>> = density_CuCrZr
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[]
NodalKernels are used for non-diffusive variables (trapped species).
[NodalKernels<<<{"href": "../../syntax/NodalKernels/index.html"}>>>]
############################## NodalKernels for W (block = 4)
[time_W]
type = TimeDerivativeNodalKernel<<<{"description": "Forms the contribution to the residual and jacobian of the time derivative term from an ODE being solved at all nodes.", "href": "../../source/nodalkernels/TimeDerivativeNodalKernel.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_W
[]
[trapping_W]
type = TrappingNodalKernel<<<{"description": "Implements a residual describing the trapping of a species in a material.", "href": "../../source/nodal_kernels/TrappingNodalKernel.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_W
temperature<<<{"description": "The temperature (K)"}>>> = temperature
alpha_t<<<{"description": "The trapping rate coefficient. This has units of 1/s (e.g. no number densities are involved)"}>>> = 2.75e11 # 1e15
N<<<{"description": "The atomic number density of the host material (1/m^3)"}>>> = 1.0e0 # = (1e0) x (${tungsten_atomic_density} #/m^3)
# Ct0 = 1.0e-4 # E.A. Hodille et al 2021 Nucl. Fusion 61 126003, trap 1
Ct0<<<{"description": "The fraction of host sites that can contribute to trapping as a function (-)"}>>> = 1.0e-4 # E.A. Hodille et al 2021 Nucl. Fusion 61 1260033, trap 2
trap_per_free<<<{"description": "An estimate for the ratio of the concentration magnitude of trapped species to free species. Setting a value for this can be helpful in producing a well-scaled matrix (-)"}>>> = 1.0e0 # 1.0e1
mobile_concentration<<<{"description": "The variable representing the mobile concentration of solute particles (1/m^3)"}>>> = 'C_mobile_W'
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[release_W]
type = ReleasingNodalKernel<<<{"description": "Implements a residual describing the release of trapped species in a material.", "href": "../../source/nodal_kernels/ReleasingNodalKernel.html"}>>>
alpha_r<<<{"description": "The release rate coefficient (1/s)"}>>> = 8.4e12 # 1.0e13
temperature<<<{"description": "The temperature (K)"}>>> = temperature
# detrapping_energy = 9863.9 # = 0.85 eV E.A. Hodille et al 2021 Nucl. Fusion 61 126003, trap 1
detrapping_energy<<<{"description": "The detrapping energy (K)"}>>> = 11604.6 # = 1.00 eV E.A. Hodille et al 2021 Nucl. Fusion 61 126003, trap 2
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_W
[]
############################## NodalKernels for Cu (block = 3)
[time_Cu]
type = TimeDerivativeNodalKernel<<<{"description": "Forms the contribution to the residual and jacobian of the time derivative term from an ODE being solved at all nodes.", "href": "../../source/nodalkernels/TimeDerivativeNodalKernel.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_Cu
[]
[trapping_Cu]
type = TrappingNodalKernel<<<{"description": "Implements a residual describing the trapping of a species in a material.", "href": "../../source/nodal_kernels/TrappingNodalKernel.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_Cu
temperature<<<{"description": "The temperature (K)"}>>> = temperature
alpha_t<<<{"description": "The trapping rate coefficient. This has units of 1/s (e.g. no number densities are involved)"}>>> = 2.75e11 # 1e15
N<<<{"description": "The atomic number density of the host material (1/m^3)"}>>> = 1.0e0 # = ${tungsten_atomic_density} #/m^3 (W lattice density)
Ct0<<<{"description": "The fraction of host sites that can contribute to trapping as a function (-)"}>>> = 5.0e-5 # R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 3
trap_per_free<<<{"description": "An estimate for the ratio of the concentration magnitude of trapped species to free species. Setting a value for this can be helpful in producing a well-scaled matrix (-)"}>>> = 1.0e0 # 1.0e1
mobile_concentration<<<{"description": "The variable representing the mobile concentration of solute particles (1/m^3)"}>>> = 'C_mobile_Cu'
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[release_Cu]
type = ReleasingNodalKernel<<<{"description": "Implements a residual describing the release of trapped species in a material.", "href": "../../source/nodal_kernels/ReleasingNodalKernel.html"}>>>
alpha_r<<<{"description": "The release rate coefficient (1/s)"}>>> = 8.4e12 # 1.0e13
temperature<<<{"description": "The temperature (K)"}>>> = temperature
detrapping_energy<<<{"description": "The detrapping energy (K)"}>>> = 5802.3 # = 0.50eV R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 3
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_Cu
[]
############################## NodalKernels for CuCrZr (block = 2)
[time_CuCrZr]
type = TimeDerivativeNodalKernel<<<{"description": "Forms the contribution to the residual and jacobian of the time derivative term from an ODE being solved at all nodes.", "href": "../../source/nodalkernels/TimeDerivativeNodalKernel.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_CuCrZr
[]
[trapping_CuCrZr]
type = TrappingNodalKernel<<<{"description": "Implements a residual describing the trapping of a species in a material.", "href": "../../source/nodal_kernels/TrappingNodalKernel.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_CuCrZr
temperature<<<{"description": "The temperature (K)"}>>> = temperature
alpha_t<<<{"description": "The trapping rate coefficient. This has units of 1/s (e.g. no number densities are involved)"}>>> = 2.75e11 # 1e15
N<<<{"description": "The atomic number density of the host material (1/m^3)"}>>> = 1.0e0 # = ${tungsten_atomic_density} #/m^3 (W lattice density)
Ct0<<<{"description": "The fraction of host sites that can contribute to trapping as a function (-)"}>>> = 5.0e-5 # R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 4
# Ct0 = 4.0e-2 # R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 5
trap_per_free<<<{"description": "An estimate for the ratio of the concentration magnitude of trapped species to free species. Setting a value for this can be helpful in producing a well-scaled matrix (-)"}>>> = 1.0e0 # 1.0e1
mobile_concentration<<<{"description": "The variable representing the mobile concentration of solute particles (1/m^3)"}>>> = 'C_mobile_CuCrZr'
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[release_CuCrZr]
type = ReleasingNodalKernel<<<{"description": "Implements a residual describing the release of trapped species in a material.", "href": "../../source/nodal_kernels/ReleasingNodalKernel.html"}>>>
alpha_r<<<{"description": "The release rate coefficient (1/s)"}>>> = 8.4e12 # 1.0e13
temperature<<<{"description": "The temperature (K)"}>>> = temperature
detrapping_energy<<<{"description": "The detrapping energy (K)"}>>> = 5802.3 # = 0.50eV R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 4
# detrapping_energy = 9631.8 # = 0.83 eV R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 5
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_CuCrZr
[]
[]
AuxVariables and AuxKernels
AuxVariables are used to track quantities that are not solved for by the PDEs, but are needed to compute materials properties, or are desirable to obtain as outputs, such as the total concentration of tritium or flux values. AuxKernels provide the expressions that define the AuxVariables
.
[AuxVariables<<<{"href": "../../syntax/AuxVariables/index.html"}>>>]
[flux_y]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[]
############################## AuxVariables for W (block = 4)
[Sc_C_mobile_W]
block = 4
[]
[Sc_C_trapped_W]
block = 4
[]
[C_total_W]
block = 4
[]
[Sc_C_total_W]
block = 4
[]
[S_empty_W]
block = 4
[]
[Sc_S_empty_W]
block = 4
[]
[S_trapped_W]
block = 4
[]
[Sc_S_trapped_W]
block = 4
[]
[S_total_W]
block = 4
[]
[Sc_S_total_W]
block = 4
[]
############################## AuxVariables for Cu (block = 3)
[Sc_C_mobile_Cu]
block = 3
[]
[Sc_C_trapped_Cu]
block = 3
[]
[C_total_Cu]
block = 3
[]
[Sc_C_total_Cu]
block = 3
[]
[S_empty_Cu]
block = 3
[]
[Sc_S_empty_Cu]
block = 3
[]
[S_trapped_Cu]
block = 3
[]
[Sc_S_trapped_Cu]
block = 3
[]
[S_total_Cu]
block = 3
[]
[Sc_S_total_Cu]
block = 3
[]
############################## AuxVariables for CuCrZr (block = 2)
[Sc_C_mobile_CuCrZr]
block = 2
[]
[Sc_C_trapped_CuCrZr]
block = 2
[]
[C_total_CuCrZr]
block = 2
[]
[Sc_C_total_CuCrZr]
block = 2
[]
[S_empty_CuCrZr]
block = 2
[]
[Sc_S_empty_CuCrZr]
block = 2
[]
[S_trapped_CuCrZr]
block = 2
[]
[Sc_S_trapped_CuCrZr]
block = 2
[]
[S_total_CuCrZr]
block = 2
[]
[Sc_S_total_CuCrZr]
block = 2
[]
[]
[AuxKernels<<<{"href": "../../syntax/AuxKernels/index.html"}>>>]
############################## AuxKernels for W (block = 4)
[Scaled_mobile_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_mobile_W
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_mobile_W
[]
[Scaled_trapped_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_trapped_W
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_trapped_W
[]
[total_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = C_total_W
type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../source/auxkernels/ParsedAux.html"}>>>
expression<<<{"description": "Parsed function expression to compute"}>>> = 'C_mobile_W + C_trapped_W'
coupled_variables<<<{"description": "Vector of coupled variable names"}>>> = 'C_mobile_W C_trapped_W'
[]
[Scaled_total_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_total_W
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_total_W
[]
[empty_sites_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_empty_W
type = EmptySitesAux<<<{"description": "Calculates the concentration of empty trapping sites.", "href": "../../source/auxkernels/EmptySitesAux.html"}>>>
N<<<{"description": "The atomic number density of the host material (1/m^3)"}>>> = '${units 1.0e0 m^-3}' # = ${tungsten_atomic_density} #/m^3 (W lattice density)
# Ct0 = ${units 1.0e-4 m^-3} # E.A. Hodille et al 2021 Nucl. Fusion 61 126003, trap 1
Ct0<<<{"description": "The fraction of host sites that can contribute to trapping as a function (-)"}>>> = '${units 1.0e-4 m^-3}' # E.A. Hodille et al 2021 Nucl. Fusion 61 1260033, trap 2
trap_per_free<<<{"description": "An estimate for the ratio of the concentration magnitude of trapped species to free species. Setting a value for this can be helpful in producing a well-scaled matrix"}>>> = 1.0e0 # 1.0e1
trapped_concentration_variables<<<{"description": "Variables representing trapped particle concentrations."}>>> = C_trapped_W
[]
[scaled_empty_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_empty_W
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_empty_W
[]
[trapped_sites_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_trapped_W
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = 1e0
source_variable<<<{"description": "The variable to be normalized"}>>> = C_trapped_W
[]
[scaled_trapped_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_trapped_W
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_trapped_W
[]
[total_sites_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_total_W
type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../source/auxkernels/ParsedAux.html"}>>>
expression<<<{"description": "Parsed function expression to compute"}>>> = 'S_trapped_W + S_empty_W'
coupled_variables<<<{"description": "Vector of coupled variable names"}>>> = 'S_trapped_W S_empty_W'
[]
[scaled_total_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_total_W
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_total_W
[]
############################## AuxKernels for Cu (block = 3)
[Scaled_mobile_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_mobile_Cu
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_mobile_Cu
[]
[Scaled_trapped_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_trapped_Cu
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_trapped_Cu
[]
[total_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = C_total_Cu
type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../source/auxkernels/ParsedAux.html"}>>>
expression<<<{"description": "Parsed function expression to compute"}>>> = 'C_mobile_Cu + C_trapped_Cu'
coupled_variables<<<{"description": "Vector of coupled variable names"}>>> = 'C_mobile_Cu C_trapped_Cu'
[]
[Scaled_total_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_total_Cu
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_total_Cu
[]
[empty_sites_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_empty_Cu
type = EmptySitesAux<<<{"description": "Calculates the concentration of empty trapping sites.", "href": "../../source/auxkernels/EmptySitesAux.html"}>>>
N<<<{"description": "The atomic number density of the host material (1/m^3)"}>>> = '${units 1.0e0 m^-3}' # = ${tungsten_atomic_density} #/m^3 (W lattice density)
Ct0<<<{"description": "The fraction of host sites that can contribute to trapping as a function (-)"}>>> = '${units 5.0e-5 m^-3}' # R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 3
trap_per_free<<<{"description": "An estimate for the ratio of the concentration magnitude of trapped species to free species. Setting a value for this can be helpful in producing a well-scaled matrix"}>>> = 1.0e0 # 1.0e1
trapped_concentration_variables<<<{"description": "Variables representing trapped particle concentrations."}>>> = C_trapped_Cu
[]
[scaled_empty_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_empty_Cu
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_empty_Cu
[]
[trapped_sites_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_trapped_Cu
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = 1e0
source_variable<<<{"description": "The variable to be normalized"}>>> = C_trapped_Cu
[]
[scaled_trapped_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_trapped_Cu
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_trapped_Cu
[]
[total_sites_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_total_Cu
type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../source/auxkernels/ParsedAux.html"}>>>
expression<<<{"description": "Parsed function expression to compute"}>>> = 'S_trapped_Cu + S_empty_Cu'
coupled_variables<<<{"description": "Vector of coupled variable names"}>>> = 'S_trapped_Cu S_empty_Cu'
[]
[scaled_total_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_total_Cu
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_total_Cu
[]
############################## AuxKernels for CuCrZr (block = 2)
[Scaled_mobile_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_mobile_CuCrZr
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_mobile_CuCrZr
[]
[Scaled_trapped_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_trapped_CuCrZr
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_trapped_CuCrZr
[]
[total_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = C_total_CuCrZr
type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../source/auxkernels/ParsedAux.html"}>>>
expression<<<{"description": "Parsed function expression to compute"}>>> = 'C_mobile_CuCrZr + C_trapped_CuCrZr'
coupled_variables<<<{"description": "Vector of coupled variable names"}>>> = 'C_mobile_CuCrZr C_trapped_CuCrZr'
[]
[Scaled_total_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_total_CuCrZr
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_total_CuCrZr
[]
[empty_sites_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_empty_CuCrZr
type = EmptySitesAux<<<{"description": "Calculates the concentration of empty trapping sites.", "href": "../../source/auxkernels/EmptySitesAux.html"}>>>
N<<<{"description": "The atomic number density of the host material (1/m^3)"}>>> = '${units 1.0e0 m^-3}' # = ${tungsten_atomic_density} #/m^3 (W lattice density)
Ct0<<<{"description": "The fraction of host sites that can contribute to trapping as a function (-)"}>>> = '${units 5.0e-5 m^-3}' # R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 4
# Ct0 = ${units 4.0e-2 m^-3} # R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 5
trap_per_free<<<{"description": "An estimate for the ratio of the concentration magnitude of trapped species to free species. Setting a value for this can be helpful in producing a well-scaled matrix"}>>> = 1.0e0 # 1.0e1
trapped_concentration_variables<<<{"description": "Variables representing trapped particle concentrations."}>>> = C_trapped_CuCrZr
[]
[scaled_empty_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_empty_CuCrZr
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_empty_CuCrZr
[]
[trapped_sites_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_trapped_CuCrZr
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = 1e0
source_variable<<<{"description": "The variable to be normalized"}>>> = C_trapped_CuCrZr
[]
[scaled_trapped_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_trapped_CuCrZr
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_trapped_CuCrZr
[]
[total_sites_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_total_CuCrZr
type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../source/auxkernels/ParsedAux.html"}>>>
expression<<<{"description": "Parsed function expression to compute"}>>> = 'S_trapped_CuCrZr + S_empty_CuCrZr'
coupled_variables<<<{"description": "Vector of coupled variable names"}>>> = 'S_trapped_CuCrZr S_empty_CuCrZr'
[]
[scaled_total_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_total_CuCrZr
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_total_CuCrZr
[]
[flux_y_W]
type = DiffusionFluxAux<<<{"description": "Compute components of flux vector for diffusion problems $(\\vec{J} = -D \\nabla C)$.", "href": "../../source/auxkernels/DiffusionFluxAux.html"}>>>
diffusivity<<<{"description": "The name of the diffusivity material property that will be used in the flux computation."}>>> = diffusivity_W
variable<<<{"description": "The name of the variable that this object applies to"}>>> = flux_y
diffusion_variable<<<{"description": "The name of the variable"}>>> = C_mobile_W
component<<<{"description": "The desired component of flux."}>>> = y
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
[]
[flux_y_Cu]
type = DiffusionFluxAux<<<{"description": "Compute components of flux vector for diffusion problems $(\\vec{J} = -D \\nabla C)$.", "href": "../../source/auxkernels/DiffusionFluxAux.html"}>>>
diffusivity<<<{"description": "The name of the diffusivity material property that will be used in the flux computation."}>>> = diffusivity_Cu
variable<<<{"description": "The name of the variable that this object applies to"}>>> = flux_y
diffusion_variable<<<{"description": "The name of the variable"}>>> = C_mobile_Cu
component<<<{"description": "The desired component of flux."}>>> = y
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
[]
[flux_y_CuCrZr]
type = DiffusionFluxAux<<<{"description": "Compute components of flux vector for diffusion problems $(\\vec{J} = -D \\nabla C)$.", "href": "../../source/auxkernels/DiffusionFluxAux.html"}>>>
diffusivity<<<{"description": "The name of the diffusivity material property that will be used in the flux computation."}>>> = diffusivity_CuCrZr
variable<<<{"description": "The name of the variable that this object applies to"}>>> = flux_y
diffusion_variable<<<{"description": "The name of the variable"}>>> = C_mobile_CuCrZr
component<<<{"description": "The desired component of flux."}>>> = y
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
[]
[]
Boundary conditions and history
For the mass conservation equation, we use FunctionNeumannBC and DirichletBC to solve the 2 term of Eq. (3). FunctionNeumannBC
is used to treat the time-dependent plasma exposure at the top (plasma-exposed) surface, and DirichletBC
is used to set the BC of the solute T atom concentration to zero at the inner CuCrZr tube (at mm). For the energy conservation equation, we use FunctionNeumannBC
and DirichletBC
to solve the 2 term of Eq. (4). FunctionNeumannBC
is used to treat the time-dependent heat flux at the top (plasma-exposed) surface, and FunctionDirichletBC is used to treat the temperature increase in the cooling tube. FunctionNeumannBC
, FunctionDirichletBC
, DirichletBC
, and FunctionDirichletBC
are MOOSE objects commonly used to represent the BCs of the variables to be solved.
We simulate a 20,000-second plasma discharge, with each 1,600-second cycle consisting of a 100-second plasma ramp-up, a 400-second steady-state plasma discharge, a 100-second plasma ramp-down, and a 1,000-second waiting phase. Up to 50 cycles are simulated to achieve the total discharge. Figure 2 shows the integrated (solute, total and trapped) tritium concentration profiles in the monoblock. It shows that implantation fluxes and temperatures vary linearly up to (from) their steady-state values from (up to) their initial values during ramp-up (ramp-down).

Figure 2: Temperature profiles (orange) and integrated tritium concentration profiles (blue) during two 1,600-second-cycle plasma discharges. This corresponds to Fig. 2 in Shimada et al. (2024).
During the steady-state plasma discharge, we set a heat flux of 10 MWm at the top of the 2D monoblock (at mm) and a cooling temperature of 552 K at the inner CuCrZr tube (at mm). We assume a 100% T plasma with a 5.0 10 ms plasma particle flux (which is half of the full 1.0 10 ms DT plasma particle flux), and only 0.1% of the incident plasma particle flux (5.0 10 ms) entered the first layer of mesh at the exposed surface ( mm) as in Shimada et al. (2024). The solute T atom concentration is set to zero at the inner CuCrZr tube (at mm).
We treated this plasma exposure by setting the flux BC of the solute T atom concentration at the exposed surface ( mm) as a simplification of the complex plasma implantation and recombination phenomena, which would require a very fine mesh and increase computational costs.
[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
[C_mob_W_top_flux]
type = FunctionNeumannBC<<<{"description": "Imposes the integrated boundary condition $\\frac{\\partial u}{\\partial n}=h(t,\\vec{x})$, where $h$ is a (possibly) time and space-dependent MOOSE Function.", "href": "../../source/bcs/FunctionNeumannBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_W
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top'
function<<<{"description": "The function."}>>> = mobile_flux_bc_func
[]
[mobile_tube]
type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../source/bcs/DirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_CuCrZr
value<<<{"description": "Value of the BC"}>>> = 1.0e-18
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = '2to1'
[]
[temp_top]
type = FunctionNeumannBC<<<{"description": "Imposes the integrated boundary condition $\\frac{\\partial u}{\\partial n}=h(t,\\vec{x})$, where $h$ is a (possibly) time and space-dependent MOOSE Function.", "href": "../../source/bcs/FunctionNeumannBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top'
function<<<{"description": "The function."}>>> = temp_flux_bc_func
[]
[temp_tube]
type = FunctionDirichletBC<<<{"description": "Imposes the essential boundary condition $u=g(t,\\vec{x})$, where $g$ is a (possibly) time and space-dependent MOOSE Function.", "href": "../../source/bcs/FunctionDirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = '2to1'
function<<<{"description": "The forcing function."}>>> = temp_inner_func
[]
[]
Functions
In this example, Functions are used to define the time-dependent tritium and heat flux at the exposed surface, as well as the coolant temperature. Functions can also be spatially-dependent. The functions are then provided to BCs, as described above.
[Functions<<<{"href": "../../syntax/Functions/index.html"}>>>]
### Maximum mobile flux of 7.90e-13 at the top surface (1.0e-4 [m])
### 1.80e23/m^2-s = (5.0e23/m^2-s) *(1-0.999) = (7.90e-13)*(${tungsten_atomic_density})/(1.0e-4) at steady state
[mobile_flux_bc_func]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../source/functions/MooseParsedFunction.html"}>>>
expression<<<{"description": "The user defined function."}>>> = 'if((t % 1600) < 100.0, 0.0 + 7.90e-13*(t % 1600)/100,
if((t % 1600) < 500.0, 7.90e-13,
if((t % 1600) < 600.0, 7.90e-13 - 7.90e-13*((t % 1600)-500)/100, 0.0)))'
[]
### Heat flux of 10MW/m^2 at steady state
[temp_flux_bc_func]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../source/functions/MooseParsedFunction.html"}>>>
expression<<<{"description": "The user defined function."}>>> = 'if((t % 1600) < 100.0, 0.0 + 1.0e7*(t % 1600)/100,
if((t % 1600) < 500.0, 1.0e7,
if((t % 1600) < 600.0, 1.0e7 - 1.0e7*((t % 1600)-500)/100, 300)))'
[]
### Maximum coolant temperature of 552K at steady state
[temp_inner_func]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../source/functions/MooseParsedFunction.html"}>>>
expression<<<{"description": "The user defined function."}>>> = 'if((t % 1600) < 100.0, 300.0 + (552-300)*(t % 1600)/100,
if((t % 1600) < 500.0, 552,
if((t % 1600) < 600.0, 552.0 - (552-300)*((t % 1600)-500)/100, 300)))'
[]
[timestep_function]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../source/functions/MooseParsedFunction.html"}>>>
expression<<<{"description": "The user defined function."}>>> = 'if((t % 1600) < 10.0, 20,
if((t % 1600) < 90.0, 40,
if((t % 1600) < 110.0, 20,
if((t % 1600) < 480.0, 40,
if((t % 1600) < 500.0, 20,
if((t % 1600) < 590.0, 4,
if((t % 1600) < 610.0, 20,
if((t % 1600) < 1500.0, 200,
if((t % 1600) < 1600.0, 40, 2)))))))))'
[]
[]
Material properties
We use the tritium mass transport properties listed in Table 2 to solve the mass conservation equation (Eq. (3)) for two variables: and (i.e., solute and trapped T atom concentrations) in three materials and the thermal properties (i.e., density, temperature-dependent specific heat, and thermal conductivity) listed in Table 3 to solve the energy conservation equation (Eq. (4)) for one variable: the temperature . The diffusivity is defined as and the solubility is defined as .
Table 2: Tritium mass transport and trapping properties used in the W, Cu, and CuCrZr layers. NOTE: Two-component solubility in W kept the maximum solubility ratio between W and Cu to 104 at low temperature. The references are provided in Shimada et al. (2024).
Material | (m/s) | (eV) | (Pa) | (eV) | Detrapping energy: (eV) | Trap density: (at.fr.) |
---|---|---|---|---|---|---|
W | 2.410 | 0.39 | 1.8710 | 1.04 | 0.85 | 1.010 |
W | 3.1410 | 0.57 | ||||
Cu | 6.610 | 0.39 | 3.1410 | 0.57 | 0.50 | 5.010 |
CuCrZr | 3.910 | 0.42 | 4.2810 | 0.39 | 0.83 | 5.010 |
Table 3: Thermal properties used in the W, Cu, and CuCrZr layers. The references are provided in Shimada et al. (2024).
Material | Density: (g m) | Specific heat: (J kg K) | Thermal conductivity: (W m K) |
---|---|---|---|
W | 19,300 | 1.1610 +7.1110 T –6.5810 T +3.2410 T –5.4510 T (293 < T (K) < 2500) | 2.4110 –2.9010 T + 2.5410 T –1.0310 T +1.5210 T (293 < T (K) < 2500) |
Cu | 8,960 | 4.2110 –6.8510 T (293 < T (K) < 873) | 3.1610 +3.1810 T –3.4910 T +1.6610 T (293 < T (K) < 873) |
CuCrZr | 8,900 | 390 | 3.8710 –1.2810 T (293 < T (K) < 927) |
[Materials<<<{"href": "../../syntax/Materials/index.html"}>>>]
############################## Materials for W (block = 4)
[diffusivity_W]
type = ADParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = diffusivity_W
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '2.4e-7*exp(-4525.8/temperature)' # H diffusivity in W
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
[solubility_W]
type = ADParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = solubility_W
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
# expression = '2.95e-5 *exp(-12069.0/temperature)' # H solubility in W = (1.87e24)/(${tungsten_atomic_density}) [#/m^3]
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '2.95e-5 *exp(-12069.0/temperature) + 4.95e-8 * exp(-6614.6/temperature)' # H solubility in W = (1.87e24)/(${tungsten_atomic_density}) [#/m^3]
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
[converter_to_regular_W]
type = MaterialADConverter<<<{"description": "Converts regular material properties to AD properties and vice versa", "href": "../../source/materials/MaterialADConverter.html"}>>>
ad_props_in<<<{"description": "The names of the AD material properties to convert to regular properties"}>>> = 'diffusivity_W'
reg_props_out<<<{"description": "The names of the output regular properties"}>>> = 'diffusivity_W_nonAD'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
[]
[heat_transfer_W]
type = GenericConstantMaterial<<<{"description": "Declares material properties based on names and values prescribed by input parameters.", "href": "../../source/materials/GenericConstantMaterial.html"}>>>
prop_names<<<{"description": "The names of the properties this material will have"}>>> = 'density_W'
prop_values<<<{"description": "The values associated with the named properties"}>>> = '19300' # [g/m^3]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
[]
[specific_heat_W]
type = ParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = specific_heat_W
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '1.16e2 + 7.11e-2 * temperature - 6.58e-5 * temperature^2 + 3.24e-8 * temperature^3 -5.45e-12 * temperature^4' # ~ 132[J/kg-K]
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
[thermal_conductivity_W]
type = ParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = thermal_conductivity_W
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
# expression = '-7.8e-9 * temperature^3 + 5.0e-5 * temperature^2 - 1.1e-1 * temperature + 1.8e2' # ~ 173.0 [ W/m-K] from R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038,
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '2.41e2 - 2.90e-1 * temperature + 2.54e-4 * temperature^2 - 1.03e-7 * temperature^3 + 1.52e-11 * temperature^4' # ~ 173.0 [ W/m-K]
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
############################## Materials for Cu (block = 3)
[diffusivity_Cu]
type = ADParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = diffusivity_Cu
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '6.60e-7*exp(-4525.8/temperature)' # H diffusivity in Cu
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
[solubility_Cu]
type = ADParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = solubility_Cu
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '4.95e-5*exp(-6614.6/temperature)' # H solubility in Cu = (3.14e24)/(${tungsten_atomic_density}) [#/m^3]
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
[converter_to_regular_Cu]
type = MaterialADConverter<<<{"description": "Converts regular material properties to AD properties and vice versa", "href": "../../source/materials/MaterialADConverter.html"}>>>
ad_props_in<<<{"description": "The names of the AD material properties to convert to regular properties"}>>> = 'diffusivity_Cu'
reg_props_out<<<{"description": "The names of the output regular properties"}>>> = 'diffusivity_Cu_nonAD'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
[]
[heat_transfer_Cu]
type = GenericConstantMaterial<<<{"description": "Declares material properties based on names and values prescribed by input parameters.", "href": "../../source/materials/GenericConstantMaterial.html"}>>>
prop_names<<<{"description": "The names of the properties this material will have"}>>> = 'density_Cu'
prop_values<<<{"description": "The values associated with the named properties"}>>> = '8960.0' # [g/m^3]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
[]
[specific_heat_Cu]
type = ParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = specific_heat_Cu
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '3.16e2 + 3.18e-1 * temperature - 3.49e-4 * temperature^2 + 1.66e-7 * temperature^3' # ~ 384 [J/kg-K]
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
[thermal_conductivity_Cu]
type = ParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = thermal_conductivity_Cu
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
# expression = '-3.9e-8 * temperature^3 + 3.8e-5 * temperature^2 - 7.9e-2 * temperature + 4.0e2' # ~ 401.0 [ W/m-K] from R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038,
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '4.21e2 - 6.85e-2 * temperature' # ~ 400.0 [ W/m-K]
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
############################## Materials for CuCrZr (block = 2)
[diffusivity_CuCrZr]
type = ADParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = diffusivity_CuCrZr
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '3.90e-7*exp(-4873.9/temperature)' # H diffusivity in CuCrZr
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
[solubility_CuCrZr]
type = ADParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = solubility_CuCrZr
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '6.75e-6*exp(-4525.8/temperature)' # H solubility in CuCrZr = (4.28e23)/(${tungsten_atomic_density}) [#/m^3]
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
[converter_to_regular_CuCrZr]
type = MaterialADConverter<<<{"description": "Converts regular material properties to AD properties and vice versa", "href": "../../source/materials/MaterialADConverter.html"}>>>
ad_props_in<<<{"description": "The names of the AD material properties to convert to regular properties"}>>> = 'diffusivity_CuCrZr'
reg_props_out<<<{"description": "The names of the output regular properties"}>>> = 'diffusivity_CuCrZr_nonAD'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
[]
[heat_transfer_CuCrZr]
type = GenericConstantMaterial<<<{"description": "Declares material properties based on names and values prescribed by input parameters.", "href": "../../source/materials/GenericConstantMaterial.html"}>>>
prop_names<<<{"description": "The names of the properties this material will have"}>>> = 'density_CuCrZr specific_heat_CuCrZr'
prop_values<<<{"description": "The values associated with the named properties"}>>> = '8900.0 390.0' # [g/m^3], [ W/m-K], [J/kg-K]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
[]
[thermal_conductivity_CuCrZr]
type = ParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = thermal_conductivity_CuCrZr
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
# expression = '5.3e-7 * temperature^3 - 6.5e-4 * temperature^2 + 2.6e-1 * temperature + 3.1e2' # ~ 320.0 [ W/m-K] from R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038,
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '3.87e2 - 1.28e-1 * temperature' # ~ 349 [ W/m-K]
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
############################## Materials for others
[interface_jump_4to3]
type = SolubilityRatioMaterial<<<{"description": "Calculates the jump in concentration across an interface.", "href": "../../source/materials/SolubilityRatioMaterial.html"}>>>
solubility_primary<<<{"description": "The material property on the primary side of the interface"}>>> = solubility_W
solubility_secondary<<<{"description": "The material property on the secondary side of the interface"}>>> = solubility_Cu
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = '4to3'
concentration_primary<<<{"description": "Primary side non-linear variable for jump computation"}>>> = C_mobile_W
concentration_secondary<<<{"description": "Secondary side non-linear variable for jump computation"}>>> = C_mobile_Cu
[]
[interface_jump_3to2]
type = SolubilityRatioMaterial<<<{"description": "Calculates the jump in concentration across an interface.", "href": "../../source/materials/SolubilityRatioMaterial.html"}>>>
solubility_primary<<<{"description": "The material property on the primary side of the interface"}>>> = solubility_Cu
solubility_secondary<<<{"description": "The material property on the secondary side of the interface"}>>> = solubility_CuCrZr
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = '3to2'
concentration_primary<<<{"description": "Primary side non-linear variable for jump computation"}>>> = C_mobile_Cu
concentration_secondary<<<{"description": "Secondary side non-linear variable for jump computation"}>>> = C_mobile_CuCrZr
[]
[]
Interface Kernels
TMAP8 assumes equilibrium between the chemical potentials of the diffusing species similar to how TMAP4 and TMAP7 treated the diffusing species across two different materials with different chemical potentials. SolubilityRatioMaterial is used to treat solute concentration differences stemming from a difference in T solubilities across the interface. The solubility ratio jump is calculated via the following: (5)
ADPenaltyInterfaceDiffusion is used to conserve the particle flux at this interface between two different solubilities. The extremely low tritium solubility in W at low temperature leads to an extremely low solute concentration in W, creating a large solute concentration difference between W and Cu. Two-component solubility in W is used to keep the maximum solubility ratio between W and Cu to 104 at low temperature to avoid the convergence issue associated with calculating two significantly different solute concentration in W and Cu.
[InterfaceKernels<<<{"href": "../../syntax/InterfaceKernels/index.html"}>>>]
[tied_4to3]
type = ADPenaltyInterfaceDiffusion<<<{"description": "A penalty-based interface condition that forcesthe continuity of variables and the flux equivalence across an interface.", "href": "../../source/interfacekernels/PenaltyInterfaceDiffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_W
neighbor_var<<<{"description": "The variable on the other side of the interface."}>>> = C_mobile_Cu
penalty<<<{"description": "The penalty that penalizes jump between primary and neighbor variables."}>>> = 0.05 # it will not converge with > 0.1, but it creates negative C_mobile _Cu with << 0.1
# jump_prop_name = solubility_ratio_4to3
jump_prop_name<<<{"description": "the name of the material property that calculates the jump."}>>> = solubility_ratio
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = '4to3'
[]
[tied_3to2]
type = ADPenaltyInterfaceDiffusion<<<{"description": "A penalty-based interface condition that forcesthe continuity of variables and the flux equivalence across an interface.", "href": "../../source/interfacekernels/PenaltyInterfaceDiffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_Cu
neighbor_var<<<{"description": "The variable on the other side of the interface."}>>> = C_mobile_CuCrZr
penalty<<<{"description": "The penalty that penalizes jump between primary and neighbor variables."}>>> = 0.05 # it will not converge with > 0.1, but it creates negative C_mobile _Cu with << 0.1
# jump_prop_name = solubility_ratio_3to2
jump_prop_name<<<{"description": "the name of the material property that calculates the jump."}>>> = solubility_ratio
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = '3to2'
[]
[]
Numerical method
We use a standard preconditioner: the “single matrix preconditioner”. The Newton method is used to model the transient tritium and thermal transport in a 2D monoblock. It is important to note that MOOSE is equipped with the built-in Message Passing Interface (MPI) protocol, as tritium and thermal transport analysis of fifty 1,600-second cycle plasma discharges in the 2D monoblock is performed in under 2 hours using a single device/computer (3.5 GHz Apple M2 Pro, 10-Core CPU/16-Core GPU) with this MPI feature.
[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
type = Transient
scheme = bdf2
solve_type = NEWTON
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
nl_rel_tol = 1e-6 # 1e-8 works for 1 cycle
nl_abs_tol = 1e-7 # 1e-11 works for 1 cycle
end_time = 8.0e4 # 50 ITER shots (3.0e4 s plasma, 2.0e4 SSP)
automatic_scaling = true
line_search = 'none'
dtmin = 1e-4
nl_max_its = 18
[TimeStepper<<<{"href": "../../syntax/Executioner/TimeStepper/index.html"}>>>]
type = IterationAdaptiveDT
dt = 20
optimal_iterations = 15
iteration_window = 1
growth_factor = 1.2
cutback_factor = 0.8
timestep_limiting_postprocessor = timestep_max_pp
[]
[]
Results
The simulation results from this example are shown in Figure 3 and Figure 4. For more results, information, and discussion about the results for this example case and their significance, the reader is referred to Shimada et al. (2024).

Figure 3: Tritium concentration profile in W (left), Cu (center), and CuCrZr (right) after ten 1,600-second cycles ( s). This corresponds to Fig. 4A in Shimada et al. (2024).

Figure 4: Tritium concentration profile in W (left), Cu (center), and CuCrZr (right) after fifty 1,600-second cycles ( s). This corresponds to Fig. 4B in Shimada et al. (2024).
gold
is a smaller version of the outputThe input file (test/tests/divertor_monoblock/divertor_monoblock.i) returns the outputs that were used in Shimada et al. (2024). However, a slightly modified version of this input is run in (test/tests/divertor_monoblock/tests) as part of TMAP8's Software Quality Assurance process: It simulates only one pulse cycle, has a coarser mesh, and outputs the results less regularly to limit the file size. As a result, the exodus file in the test gold
directory is a smaller version of the output generated when running the full input file.
Complete input file
Below is the complete input file, which can be run reliably with approximately 4 processor cores. Note that this input file has not been optimized for computational costs.
### Nomenclatures
### C_mobile_j mobile H concentration in "j" material, where j = CuCrZr, Cu, W
### C_trapped_j trapped H concentration in "j" material, where j = CuCrZr, Cu, W
### C_total_j total H concentration in "j" material, where j = CuCrZr, Cu, W
###
### S_empty_j empty site concentration in "j" material, where j = CuCrZr, Cu, W
### S_trapped_j trapped site concentration in "j" material, where j = CuCrZr, Cu, W
### S_total_j total site H concentration in "j" material, where j = CuCrZr, Cu, W
###
### F_permeation permeation flux
### F_recombination recombination flux
###
### Sc_ Scaled
### Int_ Integrated
### ScInt_ Scaled and integrated
tungsten_atomic_density = ${units 6.338e28 m^-3}
## Mesh generated by ConcentricCircleMeshGenerator and saved as "2DMonoblock.e"
[Mesh<<<{"href": "../../syntax/Mesh/index.html"}>>>]
[ccmg]
type = ConcentricCircleMeshGenerator<<<{"description": "This ConcentricCircleMeshGenerator source code is to generate concentric circle meshes.", "href": "../../source/meshgenerators/ConcentricCircleMeshGenerator.html"}>>>
num_sectors<<<{"description": "num_sectors % 2 = 0, num_sectors > 0Number of azimuthal sectors in each quadrant'num_sectors' must be an even number."}>>> = 36
rings<<<{"description": "Number of rings in each circle or in the enclosing square"}>>> = '1 30 20 110'
radii<<<{"description": "Radii of major concentric circles"}>>> = '${units 6 mm -> m} ${units 7.5 mm -> m} ${units 8.5 mm -> m}'
has_outer_square<<<{"description": "It determines if meshes for a outer square are added to concentric circle meshes."}>>> = on
pitch<<<{"description": "The enclosing square can be added to the completed concentric circle mesh.Elements are quad meshes."}>>> = ${units 28 mm -> m}
portion<<<{"description": "Control of which part of mesh is created"}>>> = left_half
preserve_volumes<<<{"description": "Volume of concentric circles can be preserved using this function."}>>> = false
smoothing_max_it<<<{"description": "Number of Laplacian smoothing iterations"}>>> = 3
[]
[ssbsg1]
type = SideSetsBetweenSubdomainsGenerator<<<{"description": "MeshGenerator that creates a sideset composed of the nodes located between two or more subdomains.", "href": "../../source/meshgenerators/SideSetsBetweenSubdomainsGenerator.html"}>>>
input<<<{"description": "The mesh we want to modify"}>>> = ccmg
primary_block<<<{"description": "The primary set of blocks for which to draw a sideset between"}>>> = '4' # W
paired_block<<<{"description": "The paired set of blocks for which to draw a sideset between"}>>> = '3' # Cu
new_boundary<<<{"description": "The list of boundary names to create on the supplied subdomain"}>>> = '4to3'
[]
[ssbsg2]
type = SideSetsBetweenSubdomainsGenerator<<<{"description": "MeshGenerator that creates a sideset composed of the nodes located between two or more subdomains.", "href": "../../source/meshgenerators/SideSetsBetweenSubdomainsGenerator.html"}>>>
input<<<{"description": "The mesh we want to modify"}>>> = ssbsg1
primary_block<<<{"description": "The primary set of blocks for which to draw a sideset between"}>>> = '3' # Cu
paired_block<<<{"description": "The paired set of blocks for which to draw a sideset between"}>>> = '4' # W
new_boundary<<<{"description": "The list of boundary names to create on the supplied subdomain"}>>> = '3to4'
[]
[ssbsg3]
type = SideSetsBetweenSubdomainsGenerator<<<{"description": "MeshGenerator that creates a sideset composed of the nodes located between two or more subdomains.", "href": "../../source/meshgenerators/SideSetsBetweenSubdomainsGenerator.html"}>>>
input<<<{"description": "The mesh we want to modify"}>>> = ssbsg2
primary_block<<<{"description": "The primary set of blocks for which to draw a sideset between"}>>> = '3' # Cu
paired_block<<<{"description": "The paired set of blocks for which to draw a sideset between"}>>> = '2' # CuCrZr
new_boundary<<<{"description": "The list of boundary names to create on the supplied subdomain"}>>> = '3to2'
[]
[ssbsg4]
type = SideSetsBetweenSubdomainsGenerator<<<{"description": "MeshGenerator that creates a sideset composed of the nodes located between two or more subdomains.", "href": "../../source/meshgenerators/SideSetsBetweenSubdomainsGenerator.html"}>>>
input<<<{"description": "The mesh we want to modify"}>>> = ssbsg3
primary_block<<<{"description": "The primary set of blocks for which to draw a sideset between"}>>> = '2' # CuCrZr
paired_block<<<{"description": "The paired set of blocks for which to draw a sideset between"}>>> = '3' # Cu
new_boundary<<<{"description": "The list of boundary names to create on the supplied subdomain"}>>> = '2to3'
[]
[ssbsg5]
type = SideSetsBetweenSubdomainsGenerator<<<{"description": "MeshGenerator that creates a sideset composed of the nodes located between two or more subdomains.", "href": "../../source/meshgenerators/SideSetsBetweenSubdomainsGenerator.html"}>>>
input<<<{"description": "The mesh we want to modify"}>>> = ssbsg4
primary_block<<<{"description": "The primary set of blocks for which to draw a sideset between"}>>> = '2' # CuCrZr
paired_block<<<{"description": "The paired set of blocks for which to draw a sideset between"}>>> = '1' # H2O
new_boundary<<<{"description": "The list of boundary names to create on the supplied subdomain"}>>> = '2to1'
[]
[bdg]
type = BlockDeletionGenerator<<<{"description": "Mesh generator which removes elements from the specified subdomains", "href": "../../source/meshgenerators/BlockDeletionGenerator.html"}>>>
input<<<{"description": "The mesh we want to modify"}>>> = ssbsg5
block<<<{"description": "The list of blocks to be deleted"}>>> = '1' # H2O
[]
[]
[Problem<<<{"href": "../../syntax/Problem/index.html"}>>>]
type = ReferenceResidualProblem
extra_tag_vectors = 'ref'
reference_vector = 'ref'
[]
[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
[temperature]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = ${units 300 K}
[]
######################### Variables for W (block = 4)
[C_mobile_W]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = ${units 1.0e-20 m^-3}
block = 4
[]
[C_trapped_W]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = ${units 1.0e-15 m^-3}
block = 4
[]
######################### Variables for Cu (block = 3)
[C_mobile_Cu]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = ${units 5.0e-17 m^-3}
block = 3
[]
[C_trapped_Cu]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = ${units 1.0e-15 m^-3}
block = 3
[]
######################### Variables for CuCrZr (block = 2)
[C_mobile_CuCrZr]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = ${units 1.0e-15 m^-3}
block = 2
[]
[C_trapped_CuCrZr]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = ${units 1.0e-15 m^-3}
block = 2
[]
[]
[AuxVariables<<<{"href": "../../syntax/AuxVariables/index.html"}>>>]
[flux_y]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[]
############################## AuxVariables for W (block = 4)
[Sc_C_mobile_W]
block = 4
[]
[Sc_C_trapped_W]
block = 4
[]
[C_total_W]
block = 4
[]
[Sc_C_total_W]
block = 4
[]
[S_empty_W]
block = 4
[]
[Sc_S_empty_W]
block = 4
[]
[S_trapped_W]
block = 4
[]
[Sc_S_trapped_W]
block = 4
[]
[S_total_W]
block = 4
[]
[Sc_S_total_W]
block = 4
[]
############################## AuxVariables for Cu (block = 3)
[Sc_C_mobile_Cu]
block = 3
[]
[Sc_C_trapped_Cu]
block = 3
[]
[C_total_Cu]
block = 3
[]
[Sc_C_total_Cu]
block = 3
[]
[S_empty_Cu]
block = 3
[]
[Sc_S_empty_Cu]
block = 3
[]
[S_trapped_Cu]
block = 3
[]
[Sc_S_trapped_Cu]
block = 3
[]
[S_total_Cu]
block = 3
[]
[Sc_S_total_Cu]
block = 3
[]
############################## AuxVariables for CuCrZr (block = 2)
[Sc_C_mobile_CuCrZr]
block = 2
[]
[Sc_C_trapped_CuCrZr]
block = 2
[]
[C_total_CuCrZr]
block = 2
[]
[Sc_C_total_CuCrZr]
block = 2
[]
[S_empty_CuCrZr]
block = 2
[]
[Sc_S_empty_CuCrZr]
block = 2
[]
[S_trapped_CuCrZr]
block = 2
[]
[Sc_S_trapped_CuCrZr]
block = 2
[]
[S_total_CuCrZr]
block = 2
[]
[Sc_S_total_CuCrZr]
block = 2
[]
[]
[Kernels<<<{"href": "../../syntax/Kernels/index.html"}>>>]
############################## Kernels for W (block = 4)
[diff_W]
type = ADMatDiffusion<<<{"description": "Diffusion equation kernel that takes an isotropic diffusivity from a material property", "href": "../../source/kernels/ADMatDiffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_W
diffusivity<<<{"description": "The diffusivity value or material property"}>>> = diffusivity_W
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[time_diff_W]
type = ADTimeDerivative<<<{"description": "The time derivative operator with the weak form of $(\\psi_i, \\frac{\\partial u_h}{\\partial t})$.", "href": "../../source/kernels/ADTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_W
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[coupled_time_W]
type = ScaledCoupledTimeDerivative<<<{"description": "Time derivative Kernel that acts on a coupled variable. Weak form: $(\\psi_i, \\frac{\\partial v_h}{\\partial t})$.", "href": "../../source/kernels/ScaledCoupledTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_W
v<<<{"description": "Coupled variable"}>>> = C_trapped_W
factor<<<{"description": "The factor by which to scale"}>>> = 1e0
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[heat_conduction_W]
type = HeatConduction<<<{"description": "Diffusive heat conduction term $-\\nabla\\cdot(k\\nabla T)$ of the thermal energy conservation equation", "href": "../../source/kernels/HeatConduction.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
diffusion_coefficient<<<{"description": "Property name of the diffusion coefficient"}>>> = thermal_conductivity_W
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[time_heat_conduction_W]
type = SpecificHeatConductionTimeDerivative<<<{"description": "Time derivative term $\\rho c_p \\frac{\\partial T}{\\partial t}$ of the heat equation with the specific heat $c_p$ and the density $\\rho$ as arguments.", "href": "../../source/kernels/SpecificHeatConductionTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
specific_heat<<<{"description": "Property name of the specific heat material property"}>>> = specific_heat_W
density<<<{"description": "Property name of the density material property"}>>> = density_W
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
############################## Kernels for Cu (block = 3)
[diff_Cu]
type = ADMatDiffusion<<<{"description": "Diffusion equation kernel that takes an isotropic diffusivity from a material property", "href": "../../source/kernels/ADMatDiffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_Cu
diffusivity<<<{"description": "The diffusivity value or material property"}>>> = diffusivity_Cu
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[time_diff_Cu]
type = ADTimeDerivative<<<{"description": "The time derivative operator with the weak form of $(\\psi_i, \\frac{\\partial u_h}{\\partial t})$.", "href": "../../source/kernels/ADTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_Cu
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[coupled_time_Cu]
type = ScaledCoupledTimeDerivative<<<{"description": "Time derivative Kernel that acts on a coupled variable. Weak form: $(\\psi_i, \\frac{\\partial v_h}{\\partial t})$.", "href": "../../source/kernels/ScaledCoupledTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_Cu
v<<<{"description": "Coupled variable"}>>> = C_trapped_Cu
factor<<<{"description": "The factor by which to scale"}>>> = 1e0
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[heat_conduction_Cu]
type = HeatConduction<<<{"description": "Diffusive heat conduction term $-\\nabla\\cdot(k\\nabla T)$ of the thermal energy conservation equation", "href": "../../source/kernels/HeatConduction.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
diffusion_coefficient<<<{"description": "Property name of the diffusion coefficient"}>>> = thermal_conductivity_Cu
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[time_heat_conduction_Cu]
type = SpecificHeatConductionTimeDerivative<<<{"description": "Time derivative term $\\rho c_p \\frac{\\partial T}{\\partial t}$ of the heat equation with the specific heat $c_p$ and the density $\\rho$ as arguments.", "href": "../../source/kernels/SpecificHeatConductionTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
specific_heat<<<{"description": "Property name of the specific heat material property"}>>> = specific_heat_Cu
density<<<{"description": "Property name of the density material property"}>>> = density_Cu
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
############################## Kernels for CuCrZr (block = 2)
[diff_CuCrZr]
type = ADMatDiffusion<<<{"description": "Diffusion equation kernel that takes an isotropic diffusivity from a material property", "href": "../../source/kernels/ADMatDiffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_CuCrZr
diffusivity<<<{"description": "The diffusivity value or material property"}>>> = diffusivity_CuCrZr
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[time_diff_CuCrZr]
type = ADTimeDerivative<<<{"description": "The time derivative operator with the weak form of $(\\psi_i, \\frac{\\partial u_h}{\\partial t})$.", "href": "../../source/kernels/ADTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_CuCrZr
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[coupled_time_CuCrZr]
type = ScaledCoupledTimeDerivative<<<{"description": "Time derivative Kernel that acts on a coupled variable. Weak form: $(\\psi_i, \\frac{\\partial v_h}{\\partial t})$.", "href": "../../source/kernels/ScaledCoupledTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_CuCrZr
v<<<{"description": "Coupled variable"}>>> = C_trapped_CuCrZr
factor<<<{"description": "The factor by which to scale"}>>> = 1e0
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[heat_conduction_CuCrZr]
type = HeatConduction<<<{"description": "Diffusive heat conduction term $-\\nabla\\cdot(k\\nabla T)$ of the thermal energy conservation equation", "href": "../../source/kernels/HeatConduction.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
diffusion_coefficient<<<{"description": "Property name of the diffusion coefficient"}>>> = thermal_conductivity_CuCrZr
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[time_heat_conduction_CuCrZr]
type = SpecificHeatConductionTimeDerivative<<<{"description": "Time derivative term $\\rho c_p \\frac{\\partial T}{\\partial t}$ of the heat equation with the specific heat $c_p$ and the density $\\rho$ as arguments.", "href": "../../source/kernels/SpecificHeatConductionTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
specific_heat<<<{"description": "Property name of the specific heat material property"}>>> = specific_heat_CuCrZr
density<<<{"description": "Property name of the density material property"}>>> = density_CuCrZr
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[]
[AuxKernels<<<{"href": "../../syntax/AuxKernels/index.html"}>>>]
############################## AuxKernels for W (block = 4)
[Scaled_mobile_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_mobile_W
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_mobile_W
[]
[Scaled_trapped_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_trapped_W
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_trapped_W
[]
[total_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = C_total_W
type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../source/auxkernels/ParsedAux.html"}>>>
expression<<<{"description": "Parsed function expression to compute"}>>> = 'C_mobile_W + C_trapped_W'
coupled_variables<<<{"description": "Vector of coupled variable names"}>>> = 'C_mobile_W C_trapped_W'
[]
[Scaled_total_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_total_W
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_total_W
[]
[empty_sites_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_empty_W
type = EmptySitesAux<<<{"description": "Calculates the concentration of empty trapping sites.", "href": "../../source/auxkernels/EmptySitesAux.html"}>>>
N<<<{"description": "The atomic number density of the host material (1/m^3)"}>>> = ${units 1.0e0 m^-3} # = ${tungsten_atomic_density} #/m^3 (W lattice density)
# Ct0 = ${units 1.0e-4 m^-3} # E.A. Hodille et al 2021 Nucl. Fusion 61 126003, trap 1
Ct0<<<{"description": "The fraction of host sites that can contribute to trapping as a function (-)"}>>> = ${units 1.0e-4 m^-3} # E.A. Hodille et al 2021 Nucl. Fusion 61 1260033, trap 2
trap_per_free<<<{"description": "An estimate for the ratio of the concentration magnitude of trapped species to free species. Setting a value for this can be helpful in producing a well-scaled matrix"}>>> = 1.0e0 # 1.0e1
trapped_concentration_variables<<<{"description": "Variables representing trapped particle concentrations."}>>> = C_trapped_W
[]
[scaled_empty_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_empty_W
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_empty_W
[]
[trapped_sites_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_trapped_W
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = 1e0
source_variable<<<{"description": "The variable to be normalized"}>>> = C_trapped_W
[]
[scaled_trapped_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_trapped_W
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_trapped_W
[]
[total_sites_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_total_W
type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../source/auxkernels/ParsedAux.html"}>>>
expression<<<{"description": "Parsed function expression to compute"}>>> = 'S_trapped_W + S_empty_W'
coupled_variables<<<{"description": "Vector of coupled variable names"}>>> = 'S_trapped_W S_empty_W'
[]
[scaled_total_W]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_total_W
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_total_W
[]
############################## AuxKernels for Cu (block = 3)
[Scaled_mobile_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_mobile_Cu
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_mobile_Cu
[]
[Scaled_trapped_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_trapped_Cu
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_trapped_Cu
[]
[total_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = C_total_Cu
type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../source/auxkernels/ParsedAux.html"}>>>
expression<<<{"description": "Parsed function expression to compute"}>>> = 'C_mobile_Cu + C_trapped_Cu'
coupled_variables<<<{"description": "Vector of coupled variable names"}>>> = 'C_mobile_Cu C_trapped_Cu'
[]
[Scaled_total_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_total_Cu
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_total_Cu
[]
[empty_sites_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_empty_Cu
type = EmptySitesAux<<<{"description": "Calculates the concentration of empty trapping sites.", "href": "../../source/auxkernels/EmptySitesAux.html"}>>>
N<<<{"description": "The atomic number density of the host material (1/m^3)"}>>> = ${units 1.0e0 m^-3} # = ${tungsten_atomic_density} #/m^3 (W lattice density)
Ct0<<<{"description": "The fraction of host sites that can contribute to trapping as a function (-)"}>>> = ${units 5.0e-5 m^-3} # R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 3
trap_per_free<<<{"description": "An estimate for the ratio of the concentration magnitude of trapped species to free species. Setting a value for this can be helpful in producing a well-scaled matrix"}>>> = 1.0e0 # 1.0e1
trapped_concentration_variables<<<{"description": "Variables representing trapped particle concentrations."}>>> = C_trapped_Cu
[]
[scaled_empty_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_empty_Cu
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_empty_Cu
[]
[trapped_sites_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_trapped_Cu
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = 1e0
source_variable<<<{"description": "The variable to be normalized"}>>> = C_trapped_Cu
[]
[scaled_trapped_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_trapped_Cu
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_trapped_Cu
[]
[total_sites_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_total_Cu
type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../source/auxkernels/ParsedAux.html"}>>>
expression<<<{"description": "Parsed function expression to compute"}>>> = 'S_trapped_Cu + S_empty_Cu'
coupled_variables<<<{"description": "Vector of coupled variable names"}>>> = 'S_trapped_Cu S_empty_Cu'
[]
[scaled_total_Cu]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_total_Cu
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_total_Cu
[]
############################## AuxKernels for CuCrZr (block = 2)
[Scaled_mobile_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_mobile_CuCrZr
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_mobile_CuCrZr
[]
[Scaled_trapped_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_trapped_CuCrZr
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_trapped_CuCrZr
[]
[total_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = C_total_CuCrZr
type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../source/auxkernels/ParsedAux.html"}>>>
expression<<<{"description": "Parsed function expression to compute"}>>> = 'C_mobile_CuCrZr + C_trapped_CuCrZr'
coupled_variables<<<{"description": "Vector of coupled variable names"}>>> = 'C_mobile_CuCrZr C_trapped_CuCrZr'
[]
[Scaled_total_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_C_total_CuCrZr
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = C_total_CuCrZr
[]
[empty_sites_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_empty_CuCrZr
type = EmptySitesAux<<<{"description": "Calculates the concentration of empty trapping sites.", "href": "../../source/auxkernels/EmptySitesAux.html"}>>>
N<<<{"description": "The atomic number density of the host material (1/m^3)"}>>> = ${units 1.0e0 m^-3} # = ${tungsten_atomic_density} #/m^3 (W lattice density)
Ct0<<<{"description": "The fraction of host sites that can contribute to trapping as a function (-)"}>>> = ${units 5.0e-5 m^-3} # R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 4
# Ct0 = ${units 4.0e-2 m^-3} # R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 5
trap_per_free<<<{"description": "An estimate for the ratio of the concentration magnitude of trapped species to free species. Setting a value for this can be helpful in producing a well-scaled matrix"}>>> = 1.0e0 # 1.0e1
trapped_concentration_variables<<<{"description": "Variables representing trapped particle concentrations."}>>> = C_trapped_CuCrZr
[]
[scaled_empty_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_empty_CuCrZr
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_empty_CuCrZr
[]
[trapped_sites_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_trapped_CuCrZr
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = 1e0
source_variable<<<{"description": "The variable to be normalized"}>>> = C_trapped_CuCrZr
[]
[scaled_trapped_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_trapped_CuCrZr
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_trapped_CuCrZr
[]
[total_sites_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = S_total_CuCrZr
type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../source/auxkernels/ParsedAux.html"}>>>
expression<<<{"description": "Parsed function expression to compute"}>>> = 'S_trapped_CuCrZr + S_empty_CuCrZr'
coupled_variables<<<{"description": "Vector of coupled variable names"}>>> = 'S_trapped_CuCrZr S_empty_CuCrZr'
[]
[scaled_total_CuCrZr]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Sc_S_total_CuCrZr
type = NormalizationAux<<<{"description": "Normalizes a variable based on a Postprocessor value.", "href": "../../source/auxkernels/NormalizationAux.html"}>>>
normal_factor<<<{"description": "The normalization factor"}>>> = ${tungsten_atomic_density}
source_variable<<<{"description": "The variable to be normalized"}>>> = S_total_CuCrZr
[]
[flux_y_W]
type = DiffusionFluxAux<<<{"description": "Compute components of flux vector for diffusion problems $(\\vec{J} = -D \\nabla C)$.", "href": "../../source/auxkernels/DiffusionFluxAux.html"}>>>
diffusivity<<<{"description": "The name of the diffusivity material property that will be used in the flux computation."}>>> = diffusivity_W
variable<<<{"description": "The name of the variable that this object applies to"}>>> = flux_y
diffusion_variable<<<{"description": "The name of the variable"}>>> = C_mobile_W
component<<<{"description": "The desired component of flux."}>>> = y
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
[]
[flux_y_Cu]
type = DiffusionFluxAux<<<{"description": "Compute components of flux vector for diffusion problems $(\\vec{J} = -D \\nabla C)$.", "href": "../../source/auxkernels/DiffusionFluxAux.html"}>>>
diffusivity<<<{"description": "The name of the diffusivity material property that will be used in the flux computation."}>>> = diffusivity_Cu
variable<<<{"description": "The name of the variable that this object applies to"}>>> = flux_y
diffusion_variable<<<{"description": "The name of the variable"}>>> = C_mobile_Cu
component<<<{"description": "The desired component of flux."}>>> = y
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
[]
[flux_y_CuCrZr]
type = DiffusionFluxAux<<<{"description": "Compute components of flux vector for diffusion problems $(\\vec{J} = -D \\nabla C)$.", "href": "../../source/auxkernels/DiffusionFluxAux.html"}>>>
diffusivity<<<{"description": "The name of the diffusivity material property that will be used in the flux computation."}>>> = diffusivity_CuCrZr
variable<<<{"description": "The name of the variable that this object applies to"}>>> = flux_y
diffusion_variable<<<{"description": "The name of the variable"}>>> = C_mobile_CuCrZr
component<<<{"description": "The desired component of flux."}>>> = y
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
[]
[]
[InterfaceKernels<<<{"href": "../../syntax/InterfaceKernels/index.html"}>>>]
[tied_4to3]
type = ADPenaltyInterfaceDiffusion<<<{"description": "A penalty-based interface condition that forcesthe continuity of variables and the flux equivalence across an interface.", "href": "../../source/interfacekernels/PenaltyInterfaceDiffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_W
neighbor_var<<<{"description": "The variable on the other side of the interface."}>>> = C_mobile_Cu
penalty<<<{"description": "The penalty that penalizes jump between primary and neighbor variables."}>>> = 0.05 # it will not converge with > 0.1, but it creates negative C_mobile _Cu with << 0.1
# jump_prop_name = solubility_ratio_4to3
jump_prop_name<<<{"description": "the name of the material property that calculates the jump."}>>> = solubility_ratio
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = '4to3'
[]
[tied_3to2]
type = ADPenaltyInterfaceDiffusion<<<{"description": "A penalty-based interface condition that forcesthe continuity of variables and the flux equivalence across an interface.", "href": "../../source/interfacekernels/PenaltyInterfaceDiffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_Cu
neighbor_var<<<{"description": "The variable on the other side of the interface."}>>> = C_mobile_CuCrZr
penalty<<<{"description": "The penalty that penalizes jump between primary and neighbor variables."}>>> = 0.05 # it will not converge with > 0.1, but it creates negative C_mobile _Cu with << 0.1
# jump_prop_name = solubility_ratio_3to2
jump_prop_name<<<{"description": "the name of the material property that calculates the jump."}>>> = solubility_ratio
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = '3to2'
[]
[]
[NodalKernels<<<{"href": "../../syntax/NodalKernels/index.html"}>>>]
############################## NodalKernels for W (block = 4)
[time_W]
type = TimeDerivativeNodalKernel<<<{"description": "Forms the contribution to the residual and jacobian of the time derivative term from an ODE being solved at all nodes.", "href": "../../source/nodalkernels/TimeDerivativeNodalKernel.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_W
[]
[trapping_W]
type = TrappingNodalKernel<<<{"description": "Implements a residual describing the trapping of a species in a material.", "href": "../../source/nodal_kernels/TrappingNodalKernel.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_W
temperature<<<{"description": "The temperature (K)"}>>> = temperature
alpha_t<<<{"description": "The trapping rate coefficient. This has units of 1/s (e.g. no number densities are involved)"}>>> = 2.75e11 # 1e15
N<<<{"description": "The atomic number density of the host material (1/m^3)"}>>> = 1.0e0 # = (1e0) x (${tungsten_atomic_density} #/m^3)
# Ct0 = 1.0e-4 # E.A. Hodille et al 2021 Nucl. Fusion 61 126003, trap 1
Ct0<<<{"description": "The fraction of host sites that can contribute to trapping as a function (-)"}>>> = 1.0e-4 # E.A. Hodille et al 2021 Nucl. Fusion 61 1260033, trap 2
trap_per_free<<<{"description": "An estimate for the ratio of the concentration magnitude of trapped species to free species. Setting a value for this can be helpful in producing a well-scaled matrix (-)"}>>> = 1.0e0 # 1.0e1
mobile_concentration<<<{"description": "The variable representing the mobile concentration of solute particles (1/m^3)"}>>> = 'C_mobile_W'
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[release_W]
type = ReleasingNodalKernel<<<{"description": "Implements a residual describing the release of trapped species in a material.", "href": "../../source/nodal_kernels/ReleasingNodalKernel.html"}>>>
alpha_r<<<{"description": "The release rate coefficient (1/s)"}>>> = 8.4e12 # 1.0e13
temperature<<<{"description": "The temperature (K)"}>>> = temperature
# detrapping_energy = 9863.9 # = 0.85 eV E.A. Hodille et al 2021 Nucl. Fusion 61 126003, trap 1
detrapping_energy<<<{"description": "The detrapping energy (K)"}>>> = 11604.6 # = 1.00 eV E.A. Hodille et al 2021 Nucl. Fusion 61 126003, trap 2
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_W
[]
############################## NodalKernels for Cu (block = 3)
[time_Cu]
type = TimeDerivativeNodalKernel<<<{"description": "Forms the contribution to the residual and jacobian of the time derivative term from an ODE being solved at all nodes.", "href": "../../source/nodalkernels/TimeDerivativeNodalKernel.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_Cu
[]
[trapping_Cu]
type = TrappingNodalKernel<<<{"description": "Implements a residual describing the trapping of a species in a material.", "href": "../../source/nodal_kernels/TrappingNodalKernel.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_Cu
temperature<<<{"description": "The temperature (K)"}>>> = temperature
alpha_t<<<{"description": "The trapping rate coefficient. This has units of 1/s (e.g. no number densities are involved)"}>>> = 2.75e11 # 1e15
N<<<{"description": "The atomic number density of the host material (1/m^3)"}>>> = 1.0e0 # = ${tungsten_atomic_density} #/m^3 (W lattice density)
Ct0<<<{"description": "The fraction of host sites that can contribute to trapping as a function (-)"}>>> = 5.0e-5 # R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 3
trap_per_free<<<{"description": "An estimate for the ratio of the concentration magnitude of trapped species to free species. Setting a value for this can be helpful in producing a well-scaled matrix (-)"}>>> = 1.0e0 # 1.0e1
mobile_concentration<<<{"description": "The variable representing the mobile concentration of solute particles (1/m^3)"}>>> = 'C_mobile_Cu'
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[release_Cu]
type = ReleasingNodalKernel<<<{"description": "Implements a residual describing the release of trapped species in a material.", "href": "../../source/nodal_kernels/ReleasingNodalKernel.html"}>>>
alpha_r<<<{"description": "The release rate coefficient (1/s)"}>>> = 8.4e12 # 1.0e13
temperature<<<{"description": "The temperature (K)"}>>> = temperature
detrapping_energy<<<{"description": "The detrapping energy (K)"}>>> = 5802.3 # = 0.50eV R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 3
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_Cu
[]
############################## NodalKernels for CuCrZr (block = 2)
[time_CuCrZr]
type = TimeDerivativeNodalKernel<<<{"description": "Forms the contribution to the residual and jacobian of the time derivative term from an ODE being solved at all nodes.", "href": "../../source/nodalkernels/TimeDerivativeNodalKernel.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_CuCrZr
[]
[trapping_CuCrZr]
type = TrappingNodalKernel<<<{"description": "Implements a residual describing the trapping of a species in a material.", "href": "../../source/nodal_kernels/TrappingNodalKernel.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_CuCrZr
temperature<<<{"description": "The temperature (K)"}>>> = temperature
alpha_t<<<{"description": "The trapping rate coefficient. This has units of 1/s (e.g. no number densities are involved)"}>>> = 2.75e11 # 1e15
N<<<{"description": "The atomic number density of the host material (1/m^3)"}>>> = 1.0e0 # = ${tungsten_atomic_density} #/m^3 (W lattice density)
Ct0<<<{"description": "The fraction of host sites that can contribute to trapping as a function (-)"}>>> = 5.0e-5 # R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 4
# Ct0 = 4.0e-2 # R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 5
trap_per_free<<<{"description": "An estimate for the ratio of the concentration magnitude of trapped species to free species. Setting a value for this can be helpful in producing a well-scaled matrix (-)"}>>> = 1.0e0 # 1.0e1
mobile_concentration<<<{"description": "The variable representing the mobile concentration of solute particles (1/m^3)"}>>> = 'C_mobile_CuCrZr'
extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = ref
[]
[release_CuCrZr]
type = ReleasingNodalKernel<<<{"description": "Implements a residual describing the release of trapped species in a material.", "href": "../../source/nodal_kernels/ReleasingNodalKernel.html"}>>>
alpha_r<<<{"description": "The release rate coefficient (1/s)"}>>> = 8.4e12 # 1.0e13
temperature<<<{"description": "The temperature (K)"}>>> = temperature
detrapping_energy<<<{"description": "The detrapping energy (K)"}>>> = 5802.3 # = 0.50eV R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 4
# detrapping_energy = 9631.8 # = 0.83 eV R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038, trap 5
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_trapped_CuCrZr
[]
[]
[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
[C_mob_W_top_flux]
type = FunctionNeumannBC<<<{"description": "Imposes the integrated boundary condition $\\frac{\\partial u}{\\partial n}=h(t,\\vec{x})$, where $h$ is a (possibly) time and space-dependent MOOSE Function.", "href": "../../source/bcs/FunctionNeumannBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_W
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top'
function<<<{"description": "The function."}>>> = mobile_flux_bc_func
[]
[mobile_tube]
type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../source/bcs/DirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = C_mobile_CuCrZr
value<<<{"description": "Value of the BC"}>>> = 1.0e-18
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = '2to1'
[]
[temp_top]
type = FunctionNeumannBC<<<{"description": "Imposes the integrated boundary condition $\\frac{\\partial u}{\\partial n}=h(t,\\vec{x})$, where $h$ is a (possibly) time and space-dependent MOOSE Function.", "href": "../../source/bcs/FunctionNeumannBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top'
function<<<{"description": "The function."}>>> = temp_flux_bc_func
[]
[temp_tube]
type = FunctionDirichletBC<<<{"description": "Imposes the essential boundary condition $u=g(t,\\vec{x})$, where $g$ is a (possibly) time and space-dependent MOOSE Function.", "href": "../../source/bcs/FunctionDirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = '2to1'
function<<<{"description": "The forcing function."}>>> = temp_inner_func
[]
[]
[Functions<<<{"href": "../../syntax/Functions/index.html"}>>>]
### Maximum mobile flux of 7.90e-13 at the top surface (1.0e-4 [m])
### 1.80e23/m^2-s = (5.0e23/m^2-s) *(1-0.999) = (7.90e-13)*(${tungsten_atomic_density})/(1.0e-4) at steady state
[mobile_flux_bc_func]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../source/functions/MooseParsedFunction.html"}>>>
expression<<<{"description": "The user defined function."}>>> = 'if((t % 1600) < 100.0, 0.0 + 7.90e-13*(t % 1600)/100,
if((t % 1600) < 500.0, 7.90e-13,
if((t % 1600) < 600.0, 7.90e-13 - 7.90e-13*((t % 1600)-500)/100, 0.0)))'
[]
### Heat flux of 10MW/m^2 at steady state
[temp_flux_bc_func]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../source/functions/MooseParsedFunction.html"}>>>
expression<<<{"description": "The user defined function."}>>> = 'if((t % 1600) < 100.0, 0.0 + 1.0e7*(t % 1600)/100,
if((t % 1600) < 500.0, 1.0e7,
if((t % 1600) < 600.0, 1.0e7 - 1.0e7*((t % 1600)-500)/100, 300)))'
[]
### Maximum coolant temperature of 552K at steady state
[temp_inner_func]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../source/functions/MooseParsedFunction.html"}>>>
expression<<<{"description": "The user defined function."}>>> = 'if((t % 1600) < 100.0, 300.0 + (552-300)*(t % 1600)/100,
if((t % 1600) < 500.0, 552,
if((t % 1600) < 600.0, 552.0 - (552-300)*((t % 1600)-500)/100, 300)))'
[]
[timestep_function]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../source/functions/MooseParsedFunction.html"}>>>
expression<<<{"description": "The user defined function."}>>> = 'if((t % 1600) < 10.0, 20,
if((t % 1600) < 90.0, 40,
if((t % 1600) < 110.0, 20,
if((t % 1600) < 480.0, 40,
if((t % 1600) < 500.0, 20,
if((t % 1600) < 590.0, 4,
if((t % 1600) < 610.0, 20,
if((t % 1600) < 1500.0, 200,
if((t % 1600) < 1600.0, 40, 2)))))))))'
[]
[]
[Materials<<<{"href": "../../syntax/Materials/index.html"}>>>]
############################## Materials for W (block = 4)
[diffusivity_W]
type = ADParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = diffusivity_W
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '2.4e-7*exp(-4525.8/temperature)' # H diffusivity in W
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
[solubility_W]
type = ADParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = solubility_W
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
# expression = '2.95e-5 *exp(-12069.0/temperature)' # H solubility in W = (1.87e24)/(${tungsten_atomic_density}) [#/m^3]
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '2.95e-5 *exp(-12069.0/temperature) + 4.95e-8 * exp(-6614.6/temperature)' # H solubility in W = (1.87e24)/(${tungsten_atomic_density}) [#/m^3]
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
[converter_to_regular_W]
type = MaterialADConverter<<<{"description": "Converts regular material properties to AD properties and vice versa", "href": "../../source/materials/MaterialADConverter.html"}>>>
ad_props_in<<<{"description": "The names of the AD material properties to convert to regular properties"}>>> = 'diffusivity_W'
reg_props_out<<<{"description": "The names of the output regular properties"}>>> = 'diffusivity_W_nonAD'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
[]
[heat_transfer_W]
type = GenericConstantMaterial<<<{"description": "Declares material properties based on names and values prescribed by input parameters.", "href": "../../source/materials/GenericConstantMaterial.html"}>>>
prop_names<<<{"description": "The names of the properties this material will have"}>>> = 'density_W'
prop_values<<<{"description": "The values associated with the named properties"}>>> = '19300' # [g/m^3]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
[]
[specific_heat_W]
type = ParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = specific_heat_W
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '1.16e2 + 7.11e-2 * temperature - 6.58e-5 * temperature^2 + 3.24e-8 * temperature^3 -5.45e-12 * temperature^4' # ~ 132[J/kg-K]
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
[thermal_conductivity_W]
type = ParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = thermal_conductivity_W
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
# expression = '-7.8e-9 * temperature^3 + 5.0e-5 * temperature^2 - 1.1e-1 * temperature + 1.8e2' # ~ 173.0 [ W/m-K] from R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038,
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '2.41e2 - 2.90e-1 * temperature + 2.54e-4 * temperature^2 - 1.03e-7 * temperature^3 + 1.52e-11 * temperature^4' # ~ 173.0 [ W/m-K]
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
############################## Materials for Cu (block = 3)
[diffusivity_Cu]
type = ADParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = diffusivity_Cu
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '6.60e-7*exp(-4525.8/temperature)' # H diffusivity in Cu
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
[solubility_Cu]
type = ADParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = solubility_Cu
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '4.95e-5*exp(-6614.6/temperature)' # H solubility in Cu = (3.14e24)/(${tungsten_atomic_density}) [#/m^3]
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
[converter_to_regular_Cu]
type = MaterialADConverter<<<{"description": "Converts regular material properties to AD properties and vice versa", "href": "../../source/materials/MaterialADConverter.html"}>>>
ad_props_in<<<{"description": "The names of the AD material properties to convert to regular properties"}>>> = 'diffusivity_Cu'
reg_props_out<<<{"description": "The names of the output regular properties"}>>> = 'diffusivity_Cu_nonAD'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
[]
[heat_transfer_Cu]
type = GenericConstantMaterial<<<{"description": "Declares material properties based on names and values prescribed by input parameters.", "href": "../../source/materials/GenericConstantMaterial.html"}>>>
prop_names<<<{"description": "The names of the properties this material will have"}>>> = 'density_Cu'
prop_values<<<{"description": "The values associated with the named properties"}>>> = '8960.0' # [g/m^3]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
[]
[specific_heat_Cu]
type = ParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = specific_heat_Cu
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '3.16e2 + 3.18e-1 * temperature - 3.49e-4 * temperature^2 + 1.66e-7 * temperature^3' # ~ 384 [J/kg-K]
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
[thermal_conductivity_Cu]
type = ParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = thermal_conductivity_Cu
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
# expression = '-3.9e-8 * temperature^3 + 3.8e-5 * temperature^2 - 7.9e-2 * temperature + 4.0e2' # ~ 401.0 [ W/m-K] from R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038,
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '4.21e2 - 6.85e-2 * temperature' # ~ 400.0 [ W/m-K]
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
############################## Materials for CuCrZr (block = 2)
[diffusivity_CuCrZr]
type = ADParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = diffusivity_CuCrZr
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '3.90e-7*exp(-4873.9/temperature)' # H diffusivity in CuCrZr
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
[solubility_CuCrZr]
type = ADParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = solubility_CuCrZr
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '6.75e-6*exp(-4525.8/temperature)' # H solubility in CuCrZr = (4.28e23)/(${tungsten_atomic_density}) [#/m^3]
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
[converter_to_regular_CuCrZr]
type = MaterialADConverter<<<{"description": "Converts regular material properties to AD properties and vice versa", "href": "../../source/materials/MaterialADConverter.html"}>>>
ad_props_in<<<{"description": "The names of the AD material properties to convert to regular properties"}>>> = 'diffusivity_CuCrZr'
reg_props_out<<<{"description": "The names of the output regular properties"}>>> = 'diffusivity_CuCrZr_nonAD'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
[]
[heat_transfer_CuCrZr]
type = GenericConstantMaterial<<<{"description": "Declares material properties based on names and values prescribed by input parameters.", "href": "../../source/materials/GenericConstantMaterial.html"}>>>
prop_names<<<{"description": "The names of the properties this material will have"}>>> = 'density_CuCrZr specific_heat_CuCrZr'
prop_values<<<{"description": "The values associated with the named properties"}>>> = '8900.0 390.0' # [g/m^3], [ W/m-K], [J/kg-K]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
[]
[thermal_conductivity_CuCrZr]
type = ParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../source/materials/ParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = thermal_conductivity_CuCrZr
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'temperature'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
# expression = '5.3e-7 * temperature^3 - 6.5e-4 * temperature^2 + 2.6e-1 * temperature + 3.1e2' # ~ 320.0 [ W/m-K] from R. Delaporte-Mathurin et al 2021 Nucl. Fusion 61 036038,
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '3.87e2 - 1.28e-1 * temperature' # ~ 349 [ W/m-K]
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = all
[]
############################## Materials for others
[interface_jump_4to3]
type = SolubilityRatioMaterial<<<{"description": "Calculates the jump in concentration across an interface.", "href": "../../source/materials/SolubilityRatioMaterial.html"}>>>
solubility_primary<<<{"description": "The material property on the primary side of the interface"}>>> = solubility_W
solubility_secondary<<<{"description": "The material property on the secondary side of the interface"}>>> = solubility_Cu
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = '4to3'
concentration_primary<<<{"description": "Primary side non-linear variable for jump computation"}>>> = C_mobile_W
concentration_secondary<<<{"description": "Secondary side non-linear variable for jump computation"}>>> = C_mobile_Cu
[]
[interface_jump_3to2]
type = SolubilityRatioMaterial<<<{"description": "Calculates the jump in concentration across an interface.", "href": "../../source/materials/SolubilityRatioMaterial.html"}>>>
solubility_primary<<<{"description": "The material property on the primary side of the interface"}>>> = solubility_Cu
solubility_secondary<<<{"description": "The material property on the secondary side of the interface"}>>> = solubility_CuCrZr
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = '3to2'
concentration_primary<<<{"description": "Primary side non-linear variable for jump computation"}>>> = C_mobile_Cu
concentration_secondary<<<{"description": "Secondary side non-linear variable for jump computation"}>>> = C_mobile_CuCrZr
[]
[]
[Postprocessors<<<{"href": "../../syntax/Postprocessors/index.html"}>>>]
############################################################ Postprocessors for W (block = 4)
[F_recombination]
type = SideDiffusiveFluxAverage<<<{"description": "Computes the integral of the diffusive flux over the specified boundary", "href": "../../source/postprocessors/SideDiffusiveFluxAverage.html"}>>>
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top'
diffusivity<<<{"description": "The name of the diffusivity material property that will be used in the flux computation. This must be provided if the variable is of finite element type"}>>> = 5.01e-24 # (3.01604928)/(6.02e23)/[gram(T)/m^2]
# diffusivity = 5.508e-19 # (1.0e3)*(1.0e3)/(6.02e23)/(3.01604928) [gram(T)/m^2]
variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = Sc_C_total_W
[]
[F_permeation]
type = SideDiffusiveFluxAverage<<<{"description": "Computes the integral of the diffusive flux over the specified boundary", "href": "../../source/postprocessors/SideDiffusiveFluxAverage.html"}>>>
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = '2to1'
diffusivity<<<{"description": "The name of the diffusivity material property that will be used in the flux computation. This must be provided if the variable is of finite element type"}>>> = 5.01e-24 # (3.01604928)/(6.02e23)/[gram(T)/m^2]
# diffusivity = 5.508e-19 # (1.0e3)*(1.0e3)/(6.02e23)/(3.01604928) [gram(T)/m^2]
variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = Sc_C_total_CuCrZr
[]
[Int_C_mobile_W]
type = ElementIntegralVariablePostprocessor<<<{"description": "Computes a volume integral of the specified variable", "href": "../../source/postprocessors/ElementIntegralVariablePostprocessor.html"}>>>
variable<<<{"description": "The name of the variable that this object operates on"}>>> = C_mobile_W
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
[]
[ScInt_C_mobile_W]
type = ScalePostprocessor<<<{"description": "Scales a post-processor by a value", "href": "../../source/postprocessors/ScalePostprocessor.html"}>>>
value<<<{"description": "The postprocessor to be scaled"}>>> = Int_C_mobile_W
scaling_factor<<<{"description": "The scaling factor"}>>> = 3.491e10 # (1.0e3)*(1.0e3)*(${tungsten_atomic_density})/(6.02e23)/(3.01604928) [gram(T)/m^2]
[]
[Int_C_trapped_W]
type = ElementIntegralVariablePostprocessor<<<{"description": "Computes a volume integral of the specified variable", "href": "../../source/postprocessors/ElementIntegralVariablePostprocessor.html"}>>>
variable<<<{"description": "The name of the variable that this object operates on"}>>> = C_trapped_W
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
[]
[ScInt_C_trapped_W]
type = ScalePostprocessor<<<{"description": "Scales a post-processor by a value", "href": "../../source/postprocessors/ScalePostprocessor.html"}>>>
value<<<{"description": "The postprocessor to be scaled"}>>> = Int_C_trapped_W
scaling_factor<<<{"description": "The scaling factor"}>>> = 3.491e10 # (1.0e3)*(1.0e3)*(${tungsten_atomic_density})/(6.02e23)/(3.01604928) [gram(T)/m^2]
[]
[Int_C_total_W]
type = ElementIntegralVariablePostprocessor<<<{"description": "Computes a volume integral of the specified variable", "href": "../../source/postprocessors/ElementIntegralVariablePostprocessor.html"}>>>
variable<<<{"description": "The name of the variable that this object operates on"}>>> = C_total_W
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 4
[]
[ScInt_C_total_W]
type = ScalePostprocessor<<<{"description": "Scales a post-processor by a value", "href": "../../source/postprocessors/ScalePostprocessor.html"}>>>
value<<<{"description": "The postprocessor to be scaled"}>>> = Int_C_total_W
scaling_factor<<<{"description": "The scaling factor"}>>> = 3.491e10 # (1.0e3)*(1.0e3)*(${tungsten_atomic_density})/(6.02e23)/(3.01604928) [gram(T)/m^2]
[]
# ############################################################ Postprocessors for Cu (block = 3)
[Int_C_mobile_Cu]
type = ElementIntegralVariablePostprocessor<<<{"description": "Computes a volume integral of the specified variable", "href": "../../source/postprocessors/ElementIntegralVariablePostprocessor.html"}>>>
variable<<<{"description": "The name of the variable that this object operates on"}>>> = C_mobile_Cu
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
[]
[ScInt_C_mobile_Cu]
type = ScalePostprocessor<<<{"description": "Scales a post-processor by a value", "href": "../../source/postprocessors/ScalePostprocessor.html"}>>>
value<<<{"description": "The postprocessor to be scaled"}>>> = Int_C_mobile_Cu
scaling_factor<<<{"description": "The scaling factor"}>>> = 3.491e10 # (1.0e3)*(1.0e3)*(${tungsten_atomic_density})/(6.02e23)/(3.01604928) [gram(T)/m^2]
[]
[Int_C_trapped_Cu]
type = ElementIntegralVariablePostprocessor<<<{"description": "Computes a volume integral of the specified variable", "href": "../../source/postprocessors/ElementIntegralVariablePostprocessor.html"}>>>
variable<<<{"description": "The name of the variable that this object operates on"}>>> = C_trapped_Cu
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
[]
[ScInt_C_trapped_Cu]
type = ScalePostprocessor<<<{"description": "Scales a post-processor by a value", "href": "../../source/postprocessors/ScalePostprocessor.html"}>>>
value<<<{"description": "The postprocessor to be scaled"}>>> = Int_C_trapped_Cu
scaling_factor<<<{"description": "The scaling factor"}>>> = 3.44e10 # (1.0e3)*(1.0e3)*(${tungsten_atomic_density})/(6.02e23)/(3.01604928) [gram(T)/m^2]
[]
[Int_C_total_Cu]
type = ElementIntegralVariablePostprocessor<<<{"description": "Computes a volume integral of the specified variable", "href": "../../source/postprocessors/ElementIntegralVariablePostprocessor.html"}>>>
variable<<<{"description": "The name of the variable that this object operates on"}>>> = C_total_Cu
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
[]
[ScInt_C_total_Cu]
type = ScalePostprocessor<<<{"description": "Scales a post-processor by a value", "href": "../../source/postprocessors/ScalePostprocessor.html"}>>>
value<<<{"description": "The postprocessor to be scaled"}>>> = Int_C_total_Cu
scaling_factor<<<{"description": "The scaling factor"}>>> = 3.491e10 # (1.0e3)*(1.0e3)*(${tungsten_atomic_density})/(6.02e23)/(3.01604928) [gram(T)/m^2]
[]
# ############################################################ Postprocessors for CuCrZr (block = 2)
[Int_C_mobile_CuCrZr]
type = ElementIntegralVariablePostprocessor<<<{"description": "Computes a volume integral of the specified variable", "href": "../../source/postprocessors/ElementIntegralVariablePostprocessor.html"}>>>
variable<<<{"description": "The name of the variable that this object operates on"}>>> = C_mobile_CuCrZr
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
[]
[ScInt_C_mobile_CuCrZr]
type = ScalePostprocessor<<<{"description": "Scales a post-processor by a value", "href": "../../source/postprocessors/ScalePostprocessor.html"}>>>
value<<<{"description": "The postprocessor to be scaled"}>>> = Int_C_mobile_CuCrZr
scaling_factor<<<{"description": "The scaling factor"}>>> = 3.491e10 # (1.0e3)*(1.0e3)*(${tungsten_atomic_density})/(6.02e23)/(3.01604928) [gram(T)/m^2]
[]
[Int_C_trapped_CuCrZr]
type = ElementIntegralVariablePostprocessor<<<{"description": "Computes a volume integral of the specified variable", "href": "../../source/postprocessors/ElementIntegralVariablePostprocessor.html"}>>>
variable<<<{"description": "The name of the variable that this object operates on"}>>> = C_trapped_CuCrZr
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
[]
[ScInt_C_trapped_CuCrZr]
type = ScalePostprocessor<<<{"description": "Scales a post-processor by a value", "href": "../../source/postprocessors/ScalePostprocessor.html"}>>>
value<<<{"description": "The postprocessor to be scaled"}>>> = Int_C_trapped_CuCrZr
scaling_factor<<<{"description": "The scaling factor"}>>> = 3.44e10 # (1.0e3)*(1.0e3)*(${tungsten_atomic_density})/(6.02e23)/(3.01604928) [gram(T)/m^2]
[]
[Int_C_total_CuCrZr]
type = ElementIntegralVariablePostprocessor<<<{"description": "Computes a volume integral of the specified variable", "href": "../../source/postprocessors/ElementIntegralVariablePostprocessor.html"}>>>
variable<<<{"description": "The name of the variable that this object operates on"}>>> = C_total_CuCrZr
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
[]
[ScInt_C_total_CuCrZr]
type = ScalePostprocessor<<<{"description": "Scales a post-processor by a value", "href": "../../source/postprocessors/ScalePostprocessor.html"}>>>
value<<<{"description": "The postprocessor to be scaled"}>>> = Int_C_total_CuCrZr
scaling_factor<<<{"description": "The scaling factor"}>>> = 3.491e10 # (1.0e3)*(1.0e3)*(${tungsten_atomic_density})/(6.02e23)/(3.01604928) [gram(T)/m^2]
[]
############################################################ Postprocessors for others
[dt]
type = TimestepSize<<<{"description": "Reports the timestep size", "href": "../../source/postprocessors/TimestepSize.html"}>>>
[]
[temperature_top]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = temperature
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0 14.0e-3 0'
[]
[temperature_tube]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = temperature
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0 6.0e-3 0'
[]
# limit timestep
[timestep_max_pp] # s
type = FunctionValuePostprocessor<<<{"description": "Computes the value of a supplied function at a single point (scalable)", "href": "../../source/postprocessors/FunctionValuePostprocessor.html"}>>>
function<<<{"description": "The function which supplies the postprocessor value."}>>> = timestep_function
[]
[]
[VectorPostprocessors<<<{"href": "../../syntax/VectorPostprocessors/index.html"}>>>]
[line]
type = LineValueSampler<<<{"description": "Samples variable(s) along a specified line", "href": "../../source/vectorpostprocessors/LineValueSampler.html"}>>>
start_point<<<{"description": "The beginning of the line"}>>> = '0 14.0e-3 0'
end_point<<<{"description": "The ending of the line"}>>> = '0 6.0e-3 0'
num_points<<<{"description": "The number of points to sample along the line"}>>> = 100
sort_by<<<{"description": "What to sort the samples by"}>>> = 'y'
variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = 'C_total_W C_total_Cu C_total_CuCrZr C_mobile_W C_mobile_Cu C_mobile_CuCrZr C_trapped_W C_trapped_Cu C_trapped_CuCrZr flux_y temperature'
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = timestep_end
[]
[]
[Preconditioning<<<{"href": "../../syntax/Preconditioning/index.html"}>>>]
[smp]
type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../../source/preconditioners/SingleMatrixPreconditioner.html"}>>>
full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
[]
[]
[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
type = Transient
scheme = bdf2
solve_type = NEWTON
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
nl_rel_tol = 1e-6 # 1e-8 works for 1 cycle
nl_abs_tol = 1e-7 # 1e-11 works for 1 cycle
end_time = 8.0e4 # 50 ITER shots (3.0e4 s plasma, 2.0e4 SSP)
automatic_scaling = true
line_search = 'none'
dtmin = 1e-4
nl_max_its = 18
[TimeStepper<<<{"href": "../../syntax/Executioner/TimeStepper/index.html"}>>>]
type = IterationAdaptiveDT
dt = 20
optimal_iterations = 15
iteration_window = 1
growth_factor = 1.2
cutback_factor = 0.8
timestep_limiting_postprocessor = timestep_max_pp
[]
[]
[Outputs<<<{"href": "../../syntax/Outputs/index.html"}>>>]
[exodus]
type = Exodus<<<{"description": "Object for output data in the Exodus format", "href": "../../source/outputs/Exodus.html"}>>>
sync_only<<<{"description": "Only export results at sync times"}>>> = false
# output at key moment in the first two cycles, and then at the end of the simulation
sync_times<<<{"description": "Times at which the output and solution is forced to occur"}>>> = '110.0 480.0 590.0 1600.0 1710.0 2080.0 2190.0 3400.0 8.0e4'
[]
csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
hide<<<{"description": "A list of the variables and postprocessors that should NOT be output to the Exodus file (may include Variables, ScalarVariables, and Postprocessor names)."}>>> = 'dt
Int_C_mobile_W Int_C_trapped_W Int_C_total_W
Int_C_mobile_Cu Int_C_trapped_Cu Int_C_total_Cu
Int_C_mobile_CuCrZr Int_C_trapped_CuCrZr Int_C_total_CuCrZr'
perf_graph<<<{"description": "Enable printing of the performance graph to the screen (Console)"}>>> = true
[]
(test/tests/divertor_monoblock/divertor_monoblock.i)References
- E. A. Hodille, R. Delaporte-Mathurin, J. Denis, M. Pecovnik, E. Bernard, Y. Ferro, R. Sakamoto, Y. Charles, J. Mougenot, A. De Backer, C. S. Becquart, S. Markelj, and C. Grisolia.
Modelling of hydrogen isotopes trapping, diffusion and permeation in divertor monoblocks under ITER-like conditions.
Nuclear Fusion, 61:126003, 10 2021.
URL: https://iopscience.iop.org/article/10.1088/1741-4326/ac2abc, doi:10.1088/1741-4326/AC2ABC.[BibTeX]
- Masashi Shimada, Pierre-Clément A. Simon, Casey T. Icenhour, and Gyanender Singh.
Toward a high-fidelity tritium transport modeling for retention and permeation experiments.
Fusion Engineering and Design, 203:114438, 6 2024.
URL: https://linkinghub.elsevier.com/retrieve/pii/S0920379624002916, doi:10.1016/J.FUSENGDES.2024.114438.[BibTeX]