discrete1.critical1d module
One-dimensional criticality drivers.
This module provides solvers for k-eigenvalue problems in one-dimensional slab and spherical geometries using power iteration methods. It includes:
Standard power iteration for multigroup criticality calculations
Data collection utilities for DJINN machine learning model training
DJINN-accelerated power iteration with surrogate model predictions
Hybrid coarse/fine mesh power iteration for improved efficiency
All drivers solve the multigroup neutron transport equation with vacuum or reflective boundary conditions to determine the critical multiplication factor (k-effective) and corresponding neutron flux distribution.
- discrete1.critical1d.collect_power_iteration(xs_total, xs_scatter, xs_fission, medium_map, delta_x, angle_x, angle_w, bc_x, filepath, geometry=1)[source]
Collect training data for DJINN models during power iteration.
Runs a standard k-eigenvalue power iteration while collecting and saving flux snapshots at each iteration step. The saved data includes scalar flux distributions, cross-section arrays, and spatial medium mapping, which can be used to train DJINN (Deep Joint Informed Neural Network) surrogate models for accelerated neutron transport calculations.
- Parameters:
xs_total (
numpy.ndarray) – Total cross sections indexed by [material, group].xs_scatter (
numpy.ndarray) – Scattering cross sections indexed by [material, group, group].xs_fission (
numpy.ndarray) – Fission cross sections indexed by [material, group].medium_map (
array_like) – Spatial material index mapping of length I (number of cells).delta_x (
array_like) – Cell widths for spatial discretization.angle_x (
array_like) – Discrete ordinates (angular quadrature points).angle_w (
array_like) – Angular quadrature weights.bc_x (
list-like) – Boundary condition indicators [left, right] (0=vacuum, 1=reflective).filepath (
str) – Directory path where training data files will be saved.geometry (
int, optional) – Geometry type (1=slab, 2=sphere). Default is 1.
- Returns:
flux (
numpy.ndarray) – Converged scalar flux distribution indexed by [cell, group].keff (
float) – Converged effective multiplication factor.
Notes
The function saves the following NumPy arrays to disk: - flux_fission_model.npy : Flux snapshots from each iteration - fission_cross_sections.npy : Fission cross section data - scatter_cross_sections.npy : Scattering cross section data - medium_map.npy : Spatial material mapping
Additionally, source iteration data is saved during the transport sweeps via mg.source_iteration_collect() for training DJINN scatter models.
- discrete1.critical1d.hybrid_power_iteration(xs_total_u, xs_scatter_u, xs_fission_u, medium_map, delta_x, angle_xu, angle_wu, bc_x, angles_c, groups_c, energy_grid, geometry=1)[source]
Run hybrid coarse-fine mesh power iteration for k-eigenvalue problems.
Implements a two-level hybrid solver that separates the solution into uncollided (fine resolution) and collided (coarse resolution) components. The fine mesh captures detailed physics while the coarse mesh handles computationally expensive collided contributions, significantly improving efficiency for problems with many energy groups and angular directions.
- Parameters:
xs_total_u (
numpy.ndarray) – Fine (uncollided) total cross sections indexed by [material, group].xs_scatter_u (
numpy.ndarray) – Fine (uncollided) scattering cross sections indexed by [material, group, group].xs_fission_u (
numpy.ndarray) – Fine (uncollided) fission cross sections indexed by [material, group].medium_map (
array_like) – Spatial material index mapping of length I (number of cells).delta_x (
array_like) – Cell widths for spatial discretization.angle_xu (
array_like) – Fine mesh discrete ordinates (angular quadrature points).angle_wu (
array_like) – Fine mesh angular quadrature weights.bc_x (
list-like) – Boundary condition indicators [left, right] (0=vacuum, 1=reflective).angles_c (
int) – Number of angular directions for coarse mesh discretization.groups_c (
int) – Number of energy groups for coarse mesh discretization.energy_grid (
tuple) – Energy grid specification (edges_g, edges_gidx_u, edges_gidx_c) defining fine-to-coarse group mapping.geometry (
int, optional) – Geometry type (1=slab, 2=sphere). Default is 1.
- Returns:
flux_u (
numpy.ndarray) – Converged fine mesh scalar flux distribution indexed by [cell, group].keff (
float) – Converged effective multiplication factor.
Notes
The hybrid method works by: 1. Computing coarse mesh collided flux with reduced angular/energy resolution 2. Projecting fission source from fine to coarse mesh 3. Combining coarse collided contribution with fine uncollided solution 4. Iterating until convergence on the fine mesh
This approach is particularly effective for deep penetration problems or systems with large spatial domains where full fine-mesh transport sweeps would be prohibitively expensive.
- discrete1.critical1d.ml_power_iteration(flux_old, xs_total, xs_scatter, xs_fission, medium_map, delta_x, angle_x, angle_w, bc_x, geometry, fission_models=[], scatter_models=[], fission_labels=None, scatter_labels=None)[source]
Power iteration with optional DJINN surrogate model acceleration.
Performs k-eigenvalue power iteration that can optionally use trained DJINN (Deep Joint Informed Neural Network) models to predict fission and/or scattering sources, significantly accelerating convergence. When no models are provided, falls back to standard power iteration.
- Parameters:
flux_old (
numpy.ndarray) – Initial scalar flux guess indexed by [cell, group].xs_total (
numpy.ndarray) – Total cross sections indexed by [material, group].xs_scatter (
numpy.ndarray) – Scattering cross sections indexed by [material, group, group].xs_fission (
numpy.ndarray) – Fission cross sections indexed by [material, group].medium_map (
array_like) – Spatial material index mapping of length I (number of cells).delta_x (
array_like) – Cell widths for spatial discretization.angle_x (
array_like) – Discrete ordinates (angular quadrature points).angle_w (
array_like) – Angular quadrature weights.bc_x (
list-like) – Boundary condition indicators [left, right] (0=vacuum, 1=reflective).geometry (
int) – Geometry type (1=slab, 2=sphere).fission_models (
list, optional) – Trained DJINN models for fission source prediction. Empty list uses traditional calculation. Default is [].scatter_models (
list, optional) – Trained DJINN models for scattering source prediction. Empty list uses traditional calculation. Default is [].fission_labels (
array_like, optional) – Material labels for fission model predictions. Default is None.scatter_labels (
array_like, optional) – Material labels for scatter model predictions. Default is None.
- Returns:
flux (
numpy.ndarray) – Converged scalar flux distribution indexed by [cell, group].keff (
float) – Converged effective multiplication factor.
Notes
Convergence criteria differ based on model usage: stricter for standard iteration (change_kk), more relaxed for DJINN-accelerated (change_pp).
Early termination occurs if flux change increases after 5 iterations, preventing divergence.
Normalization strategy changes when DJINN models are used to maintain numerical stability.
- discrete1.critical1d.power_iteration(xs_total, xs_scatter, xs_fission, medium_map, delta_x, angle_x, angle_w, bc_x, geometry=1)[source]
Run power iteration for 1D multigroup problems.
- Parameters:
xs_total (
numpy.ndarray) – Cross section arrays indexed by material.xs_scatter (
numpy.ndarray) – Cross section arrays indexed by material.xs_fission (
numpy.ndarray) – Cross section arrays indexed by material.medium_map (
array_like) – Spatial medium mapping (length I).delta_x (
array_like) – Cell widths.angle_x (
array_like) – Angular ordinates and weights.angle_w (
array_like) – Angular ordinates and weights.bc_x (
list-like) – Boundary condition indicators.geometry (
int, optional) – Geometry selector (1=slab, 2=sphere).
- Returns:
(flux, keff) converged scalar flux and multiplication factor.
- Return type:
tuple