Skip to content

Creating a new HDF5 nuclear data library for UQ#3911

Draft
Grego01-biot wants to merge 26 commits intoopenmc-dev:developfrom
Grego01-biot:covariance_storage
Draft

Creating a new HDF5 nuclear data library for UQ#3911
Grego01-biot wants to merge 26 commits intoopenmc-dev:developfrom
Grego01-biot:covariance_storage

Conversation

@Grego01-biot
Copy link
Copy Markdown
Contributor

@Grego01-biot Grego01-biot commented Apr 1, 2026

Description

The goal of this PR is to create a new HDF5 nuclear data library with covariance data stored in that would live on the OpenMC nuclear data libraries page. Another draft PR will be added shortly to openmc-dev/data repository in order to have the capability to generate the entire nuclear data library.

In order to do so, two new Python files are created and added to openmc/data and neutron.py is modified. The idea is to start by generating multigroup covariance matrices for MF33 using ERRORR module from NJOY, parsing output tape33, performing an eigenvalue decomposition (QR factorization) and storing the energy grids, the reactions and the lower triangular matrices for each reaction. The idea behind using a separate nuclear data library is that file size is less of a concern, since it is intended only for UQ purposes, and we can still reduce storage by keeping only the lower triangular matrix. This would allow the user to perform uncertainty quantification within OpenMC by doing either a sensitivity analysis with the sandwich rule or a Total Monte Carlo approach. The data structure is not fixed yet and is open to suggestions.

Here is an example of the code used to generate the new HDF5 file for Fe56.h5 using 1500 energy groups:

from __future__ import annotations
import argparse
from pathlib import Path

import h5py
import openmc.data


def uniform_lethargy_1500_grid(emin_ev: float = 1e-5, emax_ev: float = 20e6):
    import numpy as np
    return np.logspace(np.log10(emin_ev), np.log10(emax_ev), 1501)


def attach_mf33_to_openmc_h5(h5_path: Path, cov) -> None:
    with h5py.File(h5_path, "r+") as f:
        nuc_name = next(iter(f.keys()))
        cov_root = f[nuc_name].require_group("covariance")
        cov.write_mf33_group(cov_root, store_raw_covariance=True)

def main():
    p = argparse.ArgumentParser(
        description="Create one OpenMC HDF5 file from one ENDF file and attach MF=33 covariances."
    )
    p.add_argument("endf", type=Path, help="Path to ENDF file")
    p.add_argument("--njoy", required=True, help="Path to NJOY executable")
    p.add_argument("-T", "--temperature", type=float, default=293.6, help="Temperature in K")
    args = p.parse_args()

    endf_path = args.endf.resolve()
    ek = uniform_lethargy_1500_grid()

    data = openmc.data.IncidentNeutron.from_njoy(
        str(endf_path),
        temperatures=[args.temperature],
        njoy_exec=args.njoy,
    )

    out_h5 = Path.cwd() / f"{data.name}.h5"
    data.export_to_hdf5(str(out_h5), "w", libver="latest")

    from openmc.data.xs_covariance_njoy import NeutronXSCovariances
    cov = NeutronXSCovariances.from_endf(
        endf_path,
        ek,
        njoy_exec=args.njoy,
        temperature=args.temperature,
        name=data.name,
        eig_tol=1e-10,
    )

    attach_mf33_to_openmc_h5(out_h5, cov)

    print(f"Wrote: {out_h5}")

if __name__ == "__main__":
    main()

It is important to note that the code needs modular to account for the tolerance used in NJOY, the dilution cross section used for the multigroup covariance generation and the combination or not of MF32 and MF33. This capability would also allow users to generate and store their own multigroup covariance data by providing a precomputed flux spectrum as input to NJOY and selecting the desired output energy grid. Another PR will follow with the source code modifications in order to sample on the fly the neutron cross sections using the covariance data stored in the HDF5.

Fixes # (issue)

Checklist

  • I have performed a self-review of my own code
  • I have run clang-format (version 18) on any C++ source files (if applicable)
  • I have followed the style guidelines for Python source files (if applicable)
  • I have made corresponding changes to the documentation (if applicable)
  • I have added tests that prove my fix is effective or that my feature works (if applicable)

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Labels

None yet

Projects

None yet

Development

Successfully merging this pull request may close these issues.

1 participant