Lattice Tool Chains#

PWR Assembly Components#

class scarabee.reseau.FuelPin(fuel: Material, fuel_radius: float, clad: Material, clad_radius, gap: Optional[Material], gap_radius: Optional[float], num_fuel_rings: int = 1)[source]#

Represents a generic fuel pin for a PWR.

Parameters:
  • fuel (Material) – Material which describes the fuel composition, density, and temperature.

  • fuel_radius (float) – Outer radius of the fuel pellet region.

  • gap (Material, optional) – Material which describes the composition, density, and temperature of the gap between the fuel pellet and the cladding, if present.

  • gap_radius (float, optional) – Outer radius of the gap material, if present.

  • clad (Material) – Material which describes the cladding composition, density, and temperature.

  • clad_radius (float) – Outer radius of the cladding.

  • num_fuel_rings (int, default 1) – Number of rings which should be used to discretize the fuel material. Each ring will be self-shielded and depleted separately.

fuel_radius#

Outer radius of the fuel pellet region.

Type:

float

fuel_ring_materials#

Contains the Material in each fuel ring for each depletion time step.

Type:

list of list of Material

fuel_ring_flux_spectra#

Contains the average flux spectrum in each fuel ring for each depletion time step.

Type:

list of list of ndarray

fuel_dancoff_corrections#

Dancoff corrections to be used when self-shielding the fuel at each depletion time step.

Type:

list of float

gap#

Material which describes the composition, density, and temperature of the gap between the fuel pellet and the cladding, if present.

Type:

Material, optional

gap_radius#

Outer radius of the gap material, if present.

Type:

float, optional

clad#

Material which describes the cladding composition, density, and temperature.

Type:

Material

clad_radius#

Outer radius of the cladding.

Type:

float

clad_dancoff_corrections#

Dancoff corrections to be used when self-shielding the cladding at each depletion time step.

Type:

list of float

num_fuel_rings#

Number of rings which should be used to discretize the fuel material. Each ring will be self-shielded and depleted separately.

Type:

int, default 1

__init__(fuel: Material, fuel_radius: float, clad: Material, clad_radius, gap: Optional[Material], gap_radius: Optional[float], num_fuel_rings: int = 1)[source]#
load_nuclides(ndl: NDLibrary) None[source]#

Loads all the nuclides for all current materials into the data library.

Parameters:

ndl (NDLibrary) – Nuclear data library which should load the nuclides.

get_fuel_material(t: int, r: int = 0) Material[source]#

Returns the Material object for a desired fuel ring at a desired depletion time step. Ring index 0 is at the center of the pin.

Parameters:
  • t (int) – Depletion time step index.

  • r (int) – Ring index. Default is 0.

Returns:

Material defining the temperature, density, and composition for the desired fuel ring and depletion time step.

Return type:

Material

get_average_fuel_nuclide_density(t: int, nuclide: str) float[source]#

Computes the average density of a nuclide within the fuel pellet at a single depletion time step.

Parameters:
  • t (int) – Depletion time step index.

  • nuclide (str) – Name of the nuclide.

Returns:

Average density of the nuclide at depletion time step t across the fuel pellet in units of atoms per barn-cm.

Return type:

float

set_xs_for_fuel_dancoff_calculation() None[source]#

Sets the 1-group cross sections to calculate the fuel Dancoff correction.

set_xs_for_clad_dancoff_calculation(ndl: NDLibrary) None[source]#

Sets the 1-group cross sections to calculate the clad Dancoff correction.

Parameters:

ndl (NDLibrary) – Nuclear data library for obtaining potential scattering cross sections.

make_dancoff_moc_cell(moderator_xs: CrossSection, dx: float, dy: float, pintype: PinCellType, isolated: bool) SimplePinCell[source]#

Makes a simplified cell suitable for performing Dancoff correction calculations. The flat source region IDs are stored locally in the FuelPin object.

Parameters:
  • moderator_xs (CrossSection) – One group cross sections for the moderator. Total should equal absorption (i.e. no scattering) and should be equal to the macroscopic potential cross section.

  • dx (float) – Width of the cell along x.

  • dy (float) – Width of the cell along y.

  • pintype (PinCellType) – How the pin cell should be split (along x, y, or only a quadrant).

  • isolated (bool) – If True, the FSR IDs are stored for the isolated pin. Otherwise, they are stored for the full pin.

Returns:

Pin cell object for MOC Dancoff correction calculation.

Return type:

SimplifiedPinCell

populate_dancoff_fsr_indexes(isomoc: MOCDriver, fullmoc: MOCDriver) None[source]#

Obtains the flat source region indexes for all of the flat source regions used in the Dancoff correction calculations.

Parameters:
  • isomoc (MOCDriver) – MOC simulation for the isolated pin.

  • fullmoc (MOCDriver) – MOC simulation for the full geometry.

set_isolated_dancoff_fuel_sources(isomoc: MOCDriver, moderator: Material) None[source]#

Initializes the fixed sources for the isolated MOC calculation required in computing Dancoff corrections. Sources are set for a fuel Dancoff corrections calculation.

Parameters:
  • isomoc (MOCDriver) – MOC simulation for the isolated geometry.

  • moderator (Material) – Material definition for the moderator, used to obtain the potential scattering cross section.

set_isolated_dancoff_clad_sources(isomoc: MOCDriver, moderator: Material, ndl: NDLibrary) None[source]#

Initializes the fixed sources for the isolated MOC calculation required in computing Dancoff corrections. Sources are set for a clad Dancoff correction calculation.

Parameters:
  • isomoc (MOCDriver) – MOC simulation for the isolated geometry.

  • moderator (Material) – Material definition for the moderator, used to obtain the potential scattering cross section.

  • ndl (NDLibrary) – Nuclear data library for obtaining potential scattering cross sections.

set_full_dancoff_fuel_sources(fullmoc: MOCDriver, moderator: Material) None[source]#

Initializes the fixed sources for the full MOC calculation required in computing Dancoff corrections. Sources are set for a fuel Dancoff correction calculation.

Parameters:
  • fullmoc (MOCDriver) – MOC simulation for the full geometry.

  • moderator (Material) – Material definition for the moderator, used to obtain the potential scattering cross section.

set_full_dancoff_clad_sources(fullmoc: MOCDriver, moderator: Material, ndl: NDLibrary) None[source]#

Initializes the fixed sources for the full MOC calculation required in computing Dancoff corrections. Sources are set for a clad Dancoff correction calculation.

Parameters:
  • isomoc (MOCDriver) – MOC simulation for the isolated geometry.

  • moderator (Material) – Material definition for the moderator, used to obtain the potential scattering cross section.

  • ndl (NDLibrary) – Nuclear data library for obtaining potential scattering cross sections.

compute_fuel_dancoff_correction(isomoc: MOCDriver, fullmoc: MOCDriver) float[source]#

Computes the Dancoff correction for the fuel region of the fuel pin.

Parameters:
  • isomoc (MOCDriver) – MOC simulation for the isolated geometry (previously solved).

  • fullmoc (MOCDriver) – MOC simulation for the full geometry (previously solved).

Returns:

Dancoff correction for the fuel region.

Return type:

float

compute_clad_dancoff_correction(isomoc: MOCDriver, fullmoc: MOCDriver) float[source]#

Computes the Dancoff correction for the cladding region of the fuel pin.

Parameters:
  • isomoc (MOCDriver) – MOC simulation for the isolated geometry (previously solved).

  • fullmoc (MOCDriver) – MOC simulation for the full geometry (previously solved).

Returns:

Dancoff correction for the cladding region.

Return type:

float

append_fuel_dancoff_correction(C) None[source]#

Saves new Dancoff correction for the fuel that will be used for all subsequent cross section updates.

Parameters:

C (float) – New Dancoff correction.

append_clad_dancoff_correction(C) None[source]#

Saves new Dancoff correction for the cladding that will be used for all subsequent cross section updates.

Parameters:

C (float) – New Dancoff correction.

set_fuel_xs_for_depletion_step(t: int, ndl: NDLibrary) None[source]#

Constructs the CrossSection object for all fuel rings of the pin at the specified depletion step.

Parameters:
  • t (int) – Index for the depletion step.

  • ndl (NDLibrary) – Nuclear data library to use for cross sections.

set_gap_xs(ndl: NDLibrary) None[source]#

Constructs the CrossSection object for the gap between the fuel pellet and the cladding of the pin.

Parameters:

ndl (NDLibrary) – Nuclear data library to use for cross sections.

set_clad_xs_for_depletion_step(t: int, ndl: NDLibrary) None[source]#

Constructs the CrossSection object for the cladding of the pin at the specified depletion step. The depletion step only changes the Dancoff correction, not the cladding composition.

Parameters:
  • t (int) – Index for the depletion step.

  • ndl (NDLibrary) – Nuclear data library to use for cross sections.

make_moc_cell(moderator_xs: CrossSection, dx: float, dy: float, pintype: PinCellType) PinCell[source]#

Constructs the pin cell object used in for the global MOC simulation.

Parameters:
  • moderator_xs (CrossSection) – Cross sections to use for the moderator surrounding the fuel pin.

  • dx (float) – Width of the cell along x.

  • dy (float) – Width of the cell along y.

  • pintype (PinCellType) – How the pin cell should be split (along x, y, or only a quadrant).

Returns:

Pin cell suitable for the true MOC calculation.

Return type:

PinCell

populate_fsr_indexes(moc: MOCDriver) None[source]#

Obtains the flat source region indexes for all of the flat source regions used in the full MOC calculations.

Parameters:

moc (MOCDriver) – MOC simulation for the full calculations.

obtain_flux_spectra(moc: MOCDriver) None[source]#

Computes the average flux spectrum for each fuel ring from the MOC simulation. Each ring’s flux spectrum is volume averaged.

Parameters:

moc (MOCDriver) – MOC simulation for the full calculations.

compute_pin_linear_power(ndl: NDLibrary)[source]#

Computes the linear power density of the fuel pin based on the current flux spectra, in units of w / cm. Does not consider the partial pin geometry at the assembly level (i.e. a half pin in a quarter assembly).

Parameters:

ndl (NDLibrary) – Nuclear data library for the fission energy release.

Returns:

Linear power density in w / cm.

Return type:

float

normalize_flux_spectrum(f) None[source]#

Applies a multiplicative factor to the flux spectra for the fuel rings. This permits normalizing the flux to a known assembly power.

Parameters:

f (float) – Normalization factor.

predict_depletion(chain: DepletionChain, ndl: NDLibrary, dt: float, dtm1: Optional[float] = None) None[source]#

Performs the predictor in the integration of the Bateman equation. If the argument for the previous time step is not provided, CE/LI will be used. Otherwise, CE/LI is used on the first depletion step, and LE/QI is used for all subsequent time steps. The predicted material compositions are appended to the materials lists.

Parameters:
  • chain (DepletionChain) – Depletion chain to use for radioactive decay and transmutation.

  • ndl (NDLibrary) – Nuclear data library.

  • dt (float) – Durration of the time step in seconds.

  • dtm1 (float, optional) – Durration of the previous time step in seconds. Default is None.

correct_depletion(chain: DepletionChain, ndl: NDLibrary, dt: float, dtm1: Optional[float] = None) None[source]#

Performs the corrector in the integration of the Bateman equation. If the argument for the previous time step is not provided, CE/LI will be used. Otherwise, CE/LI is used on the first depletion step, and LE/QI is used for all subsequent time steps. The corrected material compositions replace the ones where were appended in the corrector step.

Parameters:
  • chain (DepletionChain) – Depletion chain to use for radioactive decay and transmutation.

  • ndl (NDLibrary) – Nuclear data library.

  • dt (float) – Durration of the time step in seconds.

  • dtm1 (float, optional) – Durration of the previous time step in seconds. Default is None.

class scarabee.reseau.GuideTube(clad: Material, inner_radius: float, outer_radius: float, fill: Optional[BurnablePoisonRod] = None)[source]#

Represents an empty guide tube for a PWR.

Parameters:
  • clad (Material) – Material which describes the cladding composition, density, and temperature.

  • inner_radius (float) – Inner radius of the guide tube.

  • outer_radius (float) – Outer radius of the guide tube.

  • fill (BurnablePoisonRod, optional) – Optional burnable poison rod which can be placed inside the guide tube. Default value is None.

clad#

Material which describes the cladding composition, density, and temperature.

Type:

Material

inner_radius#

Inner radius of the guide tube.

Type:

float

outer_radius#

Outer radius of the guide tube.

Type:

float

fill#

Optional burnable poison rod which can be placed inside the guide tube.

Type:

BurnablePoisonRod, optional

clad_dancoff_corrections#

Dancoff corrections to be used when self-shielding the cladding at each depletion time step.

Type:

list of float

empty#

True if fill is None and False otherwise.

Type:

bool

__init__(clad: Material, inner_radius: float, outer_radius: float, fill: Optional[BurnablePoisonRod] = None)[source]#
load_nuclides(ndl: NDLibrary) None[source]#

Loads all the nuclides for all current materials into the data library.

Parameters:

ndl (NDLibrary) – Nuclear data library which should load the nuclides.

set_xs_for_fuel_dancoff_calculation() None[source]#

Sets the 1-group cross sections to calculate the fuel Dancoff corrections.

set_xs_for_clad_dancoff_calculation(ndl: NDLibrary) None[source]#

Sets the 1-group cross sections to calculate the clad Dancoff corrections.

Parameters:

ndl (NDLibrary) – Nuclear data library for obtaining potential scattering cross sections.

make_dancoff_moc_cell(moderator_xs: CrossSection, dx: float, dy: float, pintype: PinCellType, isolated: bool) SimplePinCell[source]#

Makes a simplified cell suitable for performing Dancoff correction calculations. The flat source region IDs are stored locally in the GuideTube object.

Parameters:
  • moderator_xs (CrossSection) – One group cross sections for the moderator. Total should equal absorption (i.e. no scattering) and should be equal to the macroscopic potential cross section.

  • dx (float) – Width of the cell along x.

  • dy (float) – Width of the cell along y.

  • pintype (PinCellType) – How the cell should be split (along x, y, or only a quadrant).

  • isolated (bool) – If True, the FSR IDs are stored for the isolated pin. Otherwise, they are stored for the full pin.

Returns:

Pin cell object for MOC Dancoff correction calculation.

Return type:

SimplifiedPinCell

populate_dancoff_fsr_indexes(isomoc: MOCDriver, fullmoc: MOCDriver) None[source]#

Obtains the flat source region indexes for all of the flat source regions used in the Dancoff correction calculations.

Parameters:
  • isomoc (MOCDriver) – MOC simulation for the isolated pin.

  • fullmoc (MOCDriver) – MOC simulation for the full geometry.

set_isolated_dancoff_fuel_sources(isomoc: MOCDriver, moderator: Material) None[source]#

Initializes the fixed sources for the isolated MOC calculation required in computing Dancoff corrections. Sources are set for a fuel Dancoff correction calculation.

Parameters:
  • isomoc (MOCDriver) – MOC simulation for the isolated geometry.

  • moderator (Material) – Material definition for the moderator, used to obtain the potential scattering cross section.

set_isolated_dancoff_clad_sources(isomoc: MOCDriver, moderator: Material, ndl: NDLibrary) None[source]#

Initializes the fixed sources for the isolated MOC calculation required in computing Dancoff corrections. Sources are set for a clad Dancoff correction calculation.

Parameters:
  • isomoc (MOCDriver) – MOC simulation for the isolated geometry.

  • moderator (Material) – Material definition for the moderator, used to obtain the potential scattering cross section.

  • ndl (NDLibrary) – Nuclear data library for obtaining potential scattering cross sections.

set_full_dancoff_fuel_sources(fullmoc: MOCDriver, moderator: Material) None[source]#

Initializes the fixed sources for the full MOC calculation required in computing Dancoff corrections. Sources are set for a fuel Dancoff correction calculation.

Parameters:
  • fullmoc (MOCDriver) – MOC simulation for the full geometry.

  • moderator (Material) – Material definition for the moderator, used to obtain the potential scattering cross section.

set_full_dancoff_clad_sources(fullmoc: MOCDriver, moderator: Material, ndl: NDLibrary) None[source]#

Initializes the fixed sources for the full MOC calculation required in computing Dancoff corrections. Sources are set for a clad Dancoff correction calculation.

Parameters:
  • isomoc (MOCDriver) – MOC simulation for the isolated geometry.

  • moderator (Material) – Material definition for the moderator, used to obtain the potential scattering cross section.

  • ndl (NDLibrary) – Nuclear data library for obtaining potential scattering cross sections.

compute_clad_dancoff_correction(isomoc: MOCDriver, fullmoc: MOCDriver) float[source]#

Computes the Dancoff correction for the cladding region.

Parameters:
  • isomoc (MOCDriver) – MOC simulation for the isolated geometry (previously solved).

  • fullmoc (MOCDriver) – MOC simulation for the full geometry (previously solved).

Returns:

Dancoff correction for the cladding region.

Return type:

float

append_clad_dancoff_correction(C) None[source]#

Saves new Dancoff correction for the cladding that will be used for all subsequent cross section updates.

Parameters:

C (float) – New Dancoff correction.

set_clad_xs_for_depletion_step(t: int, ndl: NDLibrary) None[source]#

Constructs the CrossSection object for the cladding of the guide tube at the specified depletion step. The depletion step only changes the Dancoff correction, not the cladding composition.

Parameters:
  • t (int) – Index for the depletion step.

  • ndl (NDLibrary) – Nuclear data library to use for cross sections.

set_fill_xs_for_depletion_step(t: int, ndl: NDLibrary) None[source]#

Constructs the CrossSection objects for the fill of the guide tube at the specified depletion step. The depletion step changes the poison composition, if filled with a burnable poison rod.

Parameters:
  • t (int) – Index for the depletion step.

  • ndl (NDLibrary) – Nuclear data library to use for cross sections.

make_moc_cell(moderator_xs: CrossSection, dx: float, dy: float, pintype: PinCellType) PinCell[source]#

Constructs the pin cell object used in for the global MOC simulation.

Parameters:
  • moderator_xs (CrossSection) – Cross sections to use for the moderator surrounding the fuel pin.

  • dx (float) – Width of the cell along x.

  • dy (float) – Width of the cell along y.

  • pintype (PinCellType) – How the pin cell should be split (along x, y, or only a quadrant).

populate_fsr_indexes(moc: MOCDriver) None[source]#

Obtains the flat source region indexes for all of the flat source regions used in the full MOC calculations.

Parameters:

moc (MOCDriver) – MOC simulation for the full calculations.

obtain_flux_spectra(moc: MOCDriver) None[source]#

If the guide tube contains a burnable poison rod, the average flux spectrum in the poison is obtained from the MOC simulation.

Parameters:

moc (MOCDriver) – MOC simulation for the full calculations.

normalize_flux_spectrum(f) None[source]#

If the guide tube contains a burnable poison rod, it applies a multiplicative factor to the flux spectra for the poison. This permits normalizing the flux to a known assembly power.

Parameters:

f (float) – Normalization factor.

predict_depletion(chain: DepletionChain, ndl: NDLibrary, dt: float, dtm1: Optional[float] = None) None[source]#

Performs the predictor in the integration of the Bateman equation. If the argument for the previous time step is not provided, CE/LI will be used. Otherwise, CE/LI is used on the first depletion step, and LE/QI is used for all subsequent time steps. The predicted material compositions are appended to the materials lists.

Parameters:
  • chain (DepletionChain) – Depletion chain to use for radioactive decay and transmutation.

  • ndl (NDLibrary) – Nuclear data library.

  • dt (float) – Durration of the time step in seconds.

  • dtm1 (float, optional) – Durration of the previous time step in seconds. Default is None.

correct_depletion(chain: DepletionChain, ndl: NDLibrary, dt: float, dtm1: Optional[float] = None) None[source]#

Performs the corrector in the integration of the Bateman equation. If the argument for the previous time step is not provided, CE/LI will be used. Otherwise, CE/LI is used on the first depletion step, and LE/QI is used for all subsequent time steps. The corrected material compositions replace the ones where were appended in the corrector step.

Parameters:
  • chain (DepletionChain) – Depletion chain to use for radioactive decay and transmutation.

  • ndl (NDLibrary) – Nuclear data library.

  • dt (float) – Durration of the time step in seconds.

  • dtm1 (float, optional) – Durration of the previous time step in seconds. Default is None.

class scarabee.reseau.BurnablePoisonRod(center: Optional[Material], clad: Material, gap: Material, poison: Material, center_radius: float, inner_clad_radius: float, inner_gap_radius: float, poison_radius: float, outer_gap_radius: float, outer_clad_radius: float)[source]#

Represents a burnable poison rod in a PWR which is inserted into an empty guide tube. This can be a wet annular burnable absorber (WABA) or a borosilicate glass (BSG) burnable absorber.

Parameters:
  • center (optional Material) – Material at the center of the poison pin. If None, the center will be filled with moderator.

  • clad (Material) – Material which describes the inner and outer cladding.

  • gap (Material) – Material which describes the gap between the cladding and poison.

  • poison (Material) – Material which describes the poison.

  • center_radius (float) – Outer radius of the center of the poison rod.

  • inner_clad_radius (float) – Outer radius of the inner cladding at the center of the rod.

  • inner_gap_radius (float) – Outer radius of the inner gap between the cladding and poison.

  • poison_radius (float) – Outer radius of the poison.

  • outer_gap_radius (float) – Outer radius of the outer gap between the cladding and poison.

  • outer_clad_radius (float) – Outer radius of the outer cladding of the poison rod.

center#

Material at the center of the poison pin. If None, the center will be filled with moderator.

Type:

optional Material

clad#

Material which describes the inner and outer cladding.

Type:

Material

gap#

Material which describes the gap between the cladding and poison.

Type:

Material

poison_materials#

Contains the poison material for each depletion time step.

Type:

list of Material

center_radius#

Outer radius of the center of the poison rod.

Type:

float

inner_clad_radius#

Outer radius of the inner cladding at the center of the rod.

Type:

float

inner_gap_radius#

Outer radius of the inner gap between the cladding and poison.

Type:

float

poison_radius#

Outer radius of the poison.

Type:

float

outer_gap_radius#

Outer radius of the outer gap between the cladding and poison.

Type:

float

outer_clad_radius#

Outer radius of the outer cladding of the poison rod.

Type:

float

__init__(center: Optional[Material], clad: Material, gap: Material, poison: Material, center_radius: float, inner_clad_radius: float, inner_gap_radius: float, poison_radius: float, outer_gap_radius: float, outer_clad_radius: float)[source]#
load_nuclides(ndl: NDLibrary) None[source]#

Loads all the nuclides for all current materials into the data library.

Parameters:

ndl (NDLibrary) – Nuclear data library which should load the nuclides.

populate_dancoff_fsr_indexes(isomoc: MOCDriver, fullmoc: MOCDriver) None[source]#

Obtains the flat source region indexes for all of the flat source regions used in the Dancoff correction calculations.

Parameters:
  • isomoc (MOCDriver) – MOC simulation for the isolated pin.

  • fullmoc (MOCDriver) – MOC simulation for the full geometry.

set_xs_for_dancoff_calculation() None[source]#

Sets the cross sections for the Dancoff calculations based on the most recent poison composition.

set_isolated_dancoff_fuel_sources(isomoc: MOCDriver, moderator: Material) None[source]#

Initializes the fixed sources for the isolated MOC calculation required in computing Dancoff corrections. Sources are set for a fuel Dancoff correction calculation.

Parameters:
  • isomoc (MOCDriver) – MOC simulation for the isolated geometry.

  • moderator (Material) – Material definition for the moderator, used to obtain the potential scattering cross section.

set_isolated_dancoff_clad_sources(isomoc: MOCDriver, moderator: Material, ndl: NDLibrary) None[source]#

Initializes the fixed sources for the isolated MOC calculation required in computing Dancoff corrections. Sources are set for a clad Dancoff correction calculation.

The cladding of a burnable poison pin is no self-shielded. Therefore, this method is an alias to set_isolated_dancoff_fuel_sources.

Parameters:
  • isomoc (MOCDriver) – MOC simulation for the isolated geometry.

  • moderator (Material) – Material definition for the moderator, used to obtain the potential scattering cross section.

  • ndl (NDLibrary) – Nuclear data library for obtaining potential scattering cross sections.

set_full_dancoff_fuel_sources(fullmoc: MOCDriver, moderator: Material) None[source]#

Initializes the fixed sources for the full MOC calculation required in computing Dancoff corrections. Sources are set for a fuel Dancoff correction calculation.

Parameters:
  • fullmoc (MOCDriver) – MOC simulation for the full geometry.

  • moderator (Material) – Material definition for the moderator, used to obtain the potential scattering cross section.

set_full_dancoff_clad_sources(fullmoc: MOCDriver, moderator: Material, ndl: NDLibrary) None[source]#

Initializes the fixed sources for the full MOC calculation required in computing Dancoff corrections. Sources are set for a clad Dancoff correction calculation.

The cladding of a burnable poison pin is no self-shielded. Therefore, this method is an alias to set_isolated_dancoff_fuel_sources.

Parameters:
  • isomoc (MOCDriver) – MOC simulation for the isolated geometry.

  • moderator (Material) – Material definition for the moderator, used to obtain the potential scattering cross section.

  • ndl (NDLibrary) – Nuclear data library for obtaining potential scattering cross sections.

set_center_xs(ndl: NDLibrary) None[source]#

Constructs the CrossSection object for the material at the center of the poison rod, so long as that material is not moderator.

Parameters:

ndl (NDLibrary) – Nuclear data library to use for cross sections.

set_gap_xs(ndl: NDLibrary) None[source]#

Constructs the CrossSection object for the gap between the poison and the cladding.

Parameters:

ndl (NDLibrary) – Nuclear data library to use for cross sections.

set_clad_xs(ndl: NDLibrary) None[source]#

Constructs the CrossSection object for the cladding of the poison rod.

The cladding of a poison rod uses infinitely dilute cross sections.

Parameters:

ndl (NDLibrary) – Nuclear data library to use for cross sections.

set_poison_xs_for_depletion_step(t: int, ndl: NDLibrary) None[source]#

Constructs the CrossSection object for the poison at the specified depletion step.

Parameters:
  • t (int) – Index for the depletion step.

  • ndl (NDLibrary) – Nuclear data library to use for cross sections.

populate_fsr_indexes(moc: MOCDriver) None[source]#

Obtains the flat source region indexes for all of the flat source regions used in the full MOC calculations.

Parameters:

moc (MOCDriver) – MOC simulation for the full calculations.

obtain_flux_spectra(moc: MOCDriver) None[source]#

Computes average flux spectrum in the poison from the MOC simulation.

Parameters:

moc (MOCDriver) – MOC simulation for the full calculations.

normalize_flux_spectrum(f) None[source]#

Applies a multiplicative factor to the flux spectra for the poison. This permits normalizing the flux to a known assembly power.

Parameters:

f (float) – Normalization factor.

predict_depletion(chain: DepletionChain, ndl: NDLibrary, dt: float, dtm1: Optional[float] = None) None[source]#

Performs the predictor in the integration of the Bateman equation. If the argument for the previous time step is not provided, CE/LI will be used. Otherwise, CE/LI is used on the first depletion step, and LE/QI is used for all subsequent time steps. The predicted material compositions are appended to the materials lists.

Parameters:
  • chain (DepletionChain) – Depletion chain to use for radioactive decay and transmutation.

  • ndl (NDLibrary) – Nuclear data library.

  • dt (float) – Durration of the time step in seconds.

  • dtm1 (float, optional) – Durration of the previous time step in seconds. Default is None.

correct_depletion(chain: DepletionChain, ndl: NDLibrary, dt: float, dtm1: Optional[float] = None) None[source]#

Performs the corrector in the integration of the Bateman equation. If the argument for the previous time step is not provided, CE/LI will be used. Otherwise, CE/LI is used on the first depletion step, and LE/QI is used for all subsequent time steps. The corrected material compositions replace the ones where were appended in the corrector step.

Parameters:
  • chain (DepletionChain) – Depletion chain to use for radioactive decay and transmutation.

  • ndl (NDLibrary) – Nuclear data library.

  • dt (float) – Durration of the time step in seconds.

  • dtm1 (float, optional) – Durration of the previous time step in seconds. Default is None.

class scarabee.reseau.CriticalLeakage(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]#

Defines a critical leakage model which can be applied to an assembly.

P1 = 1#

Homogeneous P1 approximation.

B1 = 2#

Homogeneous B1 approximation.

FundamentalMode = 3#

Fundamental mode approximation.

NoLeakage = 4#

No leakage (assembly calculation is not modified).

class scarabee.reseau.Symmetry(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]#

Defines the symmetry to use when simulating a PWR fuel assembly.

Full = 1#

The fuel assembly has no symmetry.

Half = 2#

The fuel assembly has symmetry about either the x or y axis (but not both).

Quarter = 3#

The fuel assembly has symmetry about the x and y axis.

class scarabee.reseau.PWRAssembly(shape: Tuple[int, int], pitch: float, ndl: NDLibrary, cells: List[List[Union[FuelPin, GuideTube]]], moderator: dict = {}, symmetry: Symmetry = Symmetry.Quarter, independent_quadrant: bool = False, linear_power: float = 42.0, assembly_pitch: Optional[float] = None, spacer_grid_width: Optional[float] = None, spacer_grid: Optional[Material] = None, grid_sleeve_width: Optional[float] = None, grid_sleeve: Optional[Material] = None)[source]#

A PWRAssembly instance is responsible for performing all the lattice calculations necessary to produce few-group cross sections for a single PWR assembly.

Parameters:
  • shape (tuple of int) – The number of pin cells in the full assembly along x and y.

  • pitch (float) – The spacing between fuel pins.

  • ndl (NDLibrary) – Nuclear data library used for the calculation.

  • cells (list of list of FuelPin or GuideTube) – All of the cells which describe the assembly geometry. Should be consistent with the symmetry argument.

  • moderator (dict) –

    Parameters for defining the moderator. The boron concentration, temperature, and pressure will be used to call scarabee.borated_water(). Alternatively, the material key can be used to provide a user specified moderator material. Default is an empty dictionary (all default values). Acceptable keys are:

    boron-ppm:

    Moderator boron concentration in parts per million. Default is 800 if not provided. (float)

    temperature:

    Moderator temperature in Kelvin. Default is 570 if not provided. (float)

    pressure:

    Moderator pressure in MPa. Default is 15.5 if not provided. (float)

    legendre-order:

    Maximum Legendre order to load for the moderator’s anisotropic scattering data. Default is 1 if not provided. (int)

    material:

    User defined material which will be used directly. The scarabee.borated_water function will not be called. (scarabee.Material)

  • symmetry (Symmetry) – Symmetry of the fuel assembly. Default is Symmetry.Quarter.

  • independent_quadrant (bool) – If symmetry is Symmetry.Quarter and this attribute is true, the quarter assembly is treated as an independent assembly, with ADFs generated for all sides and CDFs generated for all corners. Form factors are only generated for the pins present in the quarter assembly. This improves the modeling of asymmetric burnable absorber loadings in the nodal calculation. Default value is False.

  • linear_power (float) – Linear power density of the full assembly in kW/cm. This value should not be reduced due to symmetry. Default is 42.

  • assembly_pitch (optional float) – Spacing between fuel assemblies. If None, assembly pitch is calculated from the shape and pitch. Default value is None.

  • spacer_grid_width (optional float) – Width of the spacer grid material between pin cells. Default value is None.

  • spacer_grid (optional Material) – Material defining the spacer grid between pin cells.

  • grid_sleeve_width (optional float) – Width of the grid sleeve around the assembly. Default value is None.

  • grid_sleeve (optional Material) – Material defining the grid sleeve around the assembly.

shape#

The number of pin cells in the full assembly along x and y.

Type:

tuple of int

pitch#

The spacing between fuel pins.

Type:

float

symmetry#

Symmetry of the fuel assembly.

Type:

Symmetry

independent_quadrant#

If symmetry is Symmetry.Quarter and this attribute is true, the quarter assembly is treated as an independent assembly, with ADFs generated for all sides and CDFs generated for all corners. Form factors are only generated for the pins present in the quarter assembly. This improves the modeling of asymmetric burnable absorber loadings in the nodal calculation.

Type:

bool

assembly_pitch#

Spacing between fuel assemblies.

Type:

float

fuel_volume_fraction#

Fraction of the assembly occupied by fuel.

Type:

float

moderator_volume_fraction#

Fraction of the assembly occupied by moderator.

Type:

float

linear_power#

Linear power density of the full assembly in kW/cm.

Type:

float

initial_heavy_metal_linear_mass#

Initial linear density of heavy metal in the full assembly at the beginning of life in units of kg/cm.

Type:

float

boron_ppm#

Moderator boron concentration in parts per million.

Type:

float

moderator_temp#

Moderator temperature in Kelvin.

Type:

float

moderator_pressure#

Moderator pressure in MPa.

Type:

float

moderator_legendre_order#

Maximum Legendre order to load for the moderator’s anisotropic scattering data.

Type:

int

moderator#

Material representing the assembly moderator.

Type:

Material

spacer_grid_width#

Width of the spacer grid material between pin cells.

Type:

optional float

spacer_grid#

Material defining the spacer grid between pin cells.

Type:

optional Material

grid_sleeve_width#

Width of the grid sleeve around the assembly.

Type:

optional float

grid_sleeve#

Material defining the grid sleeve around the assembly.

Type:

optional Material

dancoff_moc_track_spacing#

Spacing between tracks in the MOC calculations for determining Dancoff corrections. Default value is 0.05 cm.

Type:

float

dancoff_moc_num_angles#

Number of azimuthal angles in the MOC calculations for determining Dancoff corrections. Default value is 32.

Type:

int

dancoff_flux_tolerance#

Flux convergence tolerance for Dancoff correction calculations. Must be in range (0., 1.E-2). Default value is 1.E-5.

Type:

float

moc_track_spacing#

Spacing between tracks in the assembly MOC calculations. Default value is 0.03 cm.

Type:

float

moc_num_angles#

Number of azimuthal angles in the assembly MOC calculations. Default value is 32.

Type:

int

moc_polar_quadrature#

The polar quadrature used for the MOC calculations. Default is YamamotoTabuchi6 for isotropic calculations and Legendre8 for anisotropic calculations.

Type:

PolarQuadrature

flux_tolerance#

Flux convergence tolerance for assembly calculations. Must be in range (0., 1.E-2). Default value is 1.E-5.

Type:

float

keff_tolerance#

Keff convergence tolerance for assembly calculations. Must be in range (0., 1.E-2). Default value is 1.E-5.

Type:

float

anisotropic#

If False, isotropic scattering with the transport correction is used in the assembly calculation. Otherwise, explicit anisotropic scattering is modeled. Default value is False. If you wish to turn anisotropic scattering on, you must set this attribute before calling the solve method for the first time.

Type:

bool

cmfd#

If True, Coarse Mesh Finite Difference diffusion will be used to accelerate the assembly calculation. If True, a CMFD energy condensation scheme must be provided (see cmfd_condensation_scheme). Default value is True. If you wish to turn CMFD off/on, you must set this attribute before calling the solve method for the first time.

Type:

bool

cmfd_condensation_scheme#

Energy condensation scheme for the CMFD acceleration of the assembly calculation.

Type:

list of list of int

condensation_scheme#

Energy condensation scheme to condense from the group structure of the library to the few-groups used in the core solver.

Type:

list of list of int

prefer_moc_adf_cdf#

If True, the MOC results will be used to generate the ADFs and CDFs, even if CMFD is being used. ADFs and CDFs that are generated in this manner are typically less accurate than those generated with the CMFD results. Default value is False.

Type:

bool

leakage_corrections#

If True, the diffusion data contains few-group diffusion cross sections which DO NOT have the critical leakage model applied. Instead, a LeakageCorrections is created to update the few-group diffusion cross sections based on the in-situ leakage from the nodal diffusion solver. If False, the diffusion data does not have a LeakageCorrections instance and the few-group cross sections DO have the critical leakage model applied. Default is False.

Type:

bool

leakage_model#

Model used to determine the critical leakage flux spectrum, also known as the fundamental mode. Default method is homogeneous P1.

Type:

CriticalLeakage

depletion_exposure_steps#

1D Numpy array of assembly burn-up exposure steps, in units of MWd/kg. Default is None.

Type:

optional ndarray

depletion_time_steps#

1D Numpy array of burn-up time steps, in units of days. Default is None.

Type:

optional ndarray

corrector_transport#

If True, the corrector depletion step performs a transport calculation to obtain an updated solution for the flux based on the predictor nuclide concentrations, performing two transport calculations per time step. If False, the flux from the predictor transport calculation is used to perform the corrector step, performing only one transport calculation per time step. Default value is True.

Type:

bool

exposures#

1D Numpy array of the total assembly burn-up exposures at which material information is available, in units of MWd/kg. Default value is an empty array before solve has been called.

Type:

ndarray

times#

1D Numpy array of the total assembly burn-up times at which material information is available, in units of days. Default value is an empty array before solve has been called.

Type:

ndarray

keff#

If depletion was not performed, this is a single float with keff for the infinite assembly. If depletion was performed, this is a 1D Numpy array for the values of keff at the tabulated burn-up exposures/times. Default value is 1 before solve has been called.

Type:

float or ndarray

depletion_data#

If a single assembly calculation is being performed without depletion, this attribute will the resulting DiffusionData instance based on the few-group condensation_scheme attribute. If a depletion calculation is performed, this attribute will be a list of DiffusionData instances, with one for each burn-up point. Before solve has been called, this attribute is None.

Type:

optional DiffusionData or list of DiffusionData

__init__(shape: Tuple[int, int], pitch: float, ndl: NDLibrary, cells: List[List[Union[FuelPin, GuideTube]]], moderator: dict = {}, symmetry: Symmetry = Symmetry.Quarter, independent_quadrant: bool = False, linear_power: float = 42.0, assembly_pitch: Optional[float] = None, spacer_grid_width: Optional[float] = None, spacer_grid: Optional[Material] = None, grid_sleeve_width: Optional[float] = None, grid_sleeve: Optional[Material] = None)[source]#
get_average_fuel_nuclide_density(t: int, nuclide: str) float[source]#

Computes the average density of a nuclide within all the fuel pellets in the assembly at a single depletion time step.

Parameters:
  • t (int) – Depletion time step index.

  • nuclide (str) – Name of the nuclide.

Returns:

Average density of the nuclide at depletion time step t across all the fuel pellets in units of atoms per barn-cm.

Return type:

float

set_dancoff_moderator_xs() None[source]#

Updates the moderator cross section for all Dancoff correction calculations.

set_dancoff_spacer_grid_sleeve_xs() None[source]#

Updates the spacer grid and grid sleeve cross sections for all Dancoff correction calculations.

compute_fuel_dancoff_corrections() None[source]#

Recomputes all Dancoff corrections for the fuel regions in the problem, using the most recent material definitions. All fuel is shelf-shielded together, regardless of wether or not is is UO2 or MOX.

compute_clad_dancoff_corrections() None[source]#

Recomputes all Dancoff corrections for the fuel pin cladding regions in the problem, using the most recent material definitions. All cladding is shelf-shielded together.

apply_dancoff_corrections() None[source]#

Appends all fuel and cladding Dancoff corrections to the appropriate cell.

self_shield_and_xs_update() None[source]#

Computes a new set of Dancoff corrections for the fuel and the cladding. After, these are applied to all the cells in the problem.

plot() None[source]#

Launches the graphical geometry plotter for the assembly calculation.

set_moderator_xs() None[source]#

Updates the moderator cross section for transport calculations.

set_spacer_grid_sleeve_xs() None[source]#

Updates the spacer grid and grid sleeve cross sections for transport calculations.

recompute_all_xs() None[source]#

Computes and applies all cross sections using the most recent material information and Dancoff corrections.

recompute_all_self_shielded_xs() None[source]#

Computes and applies all fuel and caldding cross sections using the most recent material information and Dancoff corrections.

recompute_all_fuel_xs() None[source]#

Computes and applies all fuel cross sections using the most recent material information and Dancoff corrections.

recompute_all_clad_xs() None[source]#

Computes and applies all cladding cross sections using the most recent material information and Dancoff corrections.

recompute_all_gap_xs() None[source]#

Computes and applies all gap cross sections using the most recent material information.

recompute_all_guide_tube_fill_xs() None[source]#

Computes and applies all cross sections for the fill objects of guide tubes. These could be for burnable poison rods or control rods.

apply_leakage_model(scilent: bool = False) None[source]#

Applied the critical leakage model to the assembly, modifying the flux in the MOC simulation directly.

Parameters:

scilent (bool) – If True, nothing is written to the log. If False, the type of spectrum calculation, as well as kinf and the buckling are logged. Default value is False.

apply_infinite_spectrum() None[source]#

Undoes the critical flux spectrum adjustment to the MOCDriver. This permits subsequent transport calcualtions to converge much faster.

obtain_flux_spectra() None[source]#

Computes the average flux spectrum for material regions whcih are depleted. This includes each fuel ring in fuel pins and the poison in burnable poison rods.

normalize_flux_to_power() None[source]#

Normalizes the flux spectra based on the specified linear power density for the assembly. It assumes that all power comes from fission.

solve() None[source]#

Solves the assembly problem. If depletion_exposure_steps or depletion_time_steps are None, then a single k-eigenvalue problem will be run. Otherwise, the depletion calculation will be performed.

class scarabee.reseau.Reflector(fuel: CrossSection, moderator: CrossSection, assembly_width: float, gap_width: float, baffle_width: float, baffle: Material, ndl: NDLibrary)[source]#

A Reflector instance is responsible for performing transport calculations necessary to produce few-group cross sections for the reflector of an LWR. The core baffle cross sections are self-shielded as an infinite slab, using the Roman two-term approximation. The calculation is performed using a 1D Sn simulation, in the group structure of the nuclear data library. This removes the need to obtain a fine-group spectrum that would be used to condense to an intermediate group structure.

Parameters:
  • fuel (CrossSection) – Homogenized cross section which is representative of a pin cell. This is typically obtained from a previous lattice calcualtion.

  • moderator (CrossSection) – Material cross sections for the moderator at desired temperature and density.

  • assembly_width (float) – Width of a single fuel assembly (and the reflector to be modeled).

  • gap_width (float) – Width of the moderator gap between the assembly and the core baffle.

  • baffle_width (float) – Width of the core baffle.

  • baffle (Material) – Material for the core baffle at desired temperature and density.

  • ndl (NDLibrary) – Nuclear data library for constructing the baffle cross sections.

condensation_scheme#

Defines how the energy groups will be condensed to the few-group structure used in nodal calculations.

Type:

list of pairs of ints

fuel#

Cross sections for a homogenized fuel assembly.

Type:

CrossSection

moderator#

Cross sections for the moderator.

Type:

CrossSection

assembly_width#

Width of fuel assembly and reflector.

Type:

float

gap_width#

Width of the moderator gap between a fuel assembly and the core baffle.

Type:

float

baffle_width#

Width of the core baffle.

Type:

float

baffle#

Self-shielded cross sections for the core baffle.

Type:

CrossSection

nangles#

Number of angles used in the Sn solver. Default is 16. Must be one of: 2, 4, 6, 8, 10, 12, 14, 16, 32, 64, 128.

Type:

int

anisotropic#

If True, the reflector calculation is performed with explicit anisotropic scattering. Otherwise, the transport correction with isotropic scattering is used. Default value is False.

Type:

bool

keff_tolerance#

Convergence criteria for keff. Default is 1.E-5.

Type:

float

flux_tolerance#

Convergence criteria for the flux. Default is 1.E-5.

Type:

float

diffusion_xs#

The few-group diffusion group constants for the reflector.

Type:

DiffusionCrossSection

adf#

The assembly discontinuity factors.

Type:

ndarray

diffusion_data#

The few-group diffusion cross sections and ADFs for the reflector.

Type:

DiffusionData

__init__(fuel: CrossSection, moderator: CrossSection, assembly_width: float, gap_width: float, baffle_width: float, baffle: Material, ndl: NDLibrary)[source]#
solve() None[source]#

Runs a 1D problem to generate few group cross sections for the reflector, with the core baffle.

Equivalence Theory#

class scarabee.reseau.NodalFlux1D(x_min: float, x_max: float, keff: float, xs: DiffusionCrossSection, avg_flx: ndarray, j_neg: ndarray, j_pos: ndarray)[source]#
__init__(x_min: float, x_max: float, keff: float, xs: DiffusionCrossSection, avg_flx: ndarray, j_neg: ndarray, j_pos: ndarray)[source]#

Performs the 1D nodal diffusion calculation using the Nodal Expansion Method (NEM) with reference currents and average fluxes. This represents the transverse integrated homogeneous nodal flux solution. The primary use of this class is for computing the homogeneous flux that is required to define discontinuity factors for generalized equivalence theory.

Parameters:
  • x_min (float) – Position of the negative surface of the node.

  • x_max (float) – Position of the positive surface of the node.

  • keff (float) – Multiplication factor from reference Sn calculation.

  • xs (DiffusionCrossSection) – Diffusion cross sections homogenized for the node.

  • avg_flx (np.ndarray) – Average flux in each group within the node from reference Sn calculation.

  • j_neg (np.ndarray) – Reference net current in each group on the negative boundary from the reference calculation.

  • j_pos (np.ndarray) – Reference net current in each group on the positive boundary from the reference calculation.

__call__(x: float | numpy.ndarray, g: int | None = None) float | numpy.ndarray[source]#

Evaluates the nodal flux at position x in group g.

Parameters:
  • x (float | np.ndarray) – Position or array of positions for the flux evaluation.

  • g (int | None) – Energy group for the flux evaluation. If none, array for all groups is returned.

Returns:

Nodal flux in group g, or in all groups

Return type:

float | np.ndarray

pos_surf_flux(g: Optional[int] = None) Union[ndarray, float][source]#

Evaluates the nodal flux at the positive surface in group g. If g is None, an array with the surface flux in all groups is returned.

Parameters:

g (int, optional) – Energy group for the flux evaluation. Default is None.

Returns:

Nodal flux on the positive surface

Return type:

np.ndarray or float

neg_surf_flux(g: Optional[int] = None) Union[ndarray, float][source]#

Evaluates the nodal flux at the negative surface in group g. If g is None, an array with the surface flux in all groups is returned.

Parameters:

g (int, optional) – Energy group for the flux evaluation. Default is None.

Returns:

Nodal flux on the negative surface

Return type:

np.ndarray or float

class scarabee.reseau.NodalFlux2D(dx: float, dy: float, keff: float, xs: DiffusionCrossSection, avg_flx: ndarray, j_x_neg: ndarray, j_x_pos: ndarray, j_y_neg: ndarray, j_y_pos: ndarray)[source]#
__init__(dx: float, dy: float, keff: float, xs: DiffusionCrossSection, avg_flx: ndarray, j_x_neg: ndarray, j_x_pos: ndarray, j_y_neg: ndarray, j_y_pos: ndarray)[source]#

Performs the 2D nodal diffusion calculation using the Nodal Expansion Method (NEM) with reference currents and average fluxes. This represents the homogeneous nodal flux solution. The primary use of this class is for computing the homogeneous flux that is required to define corner discontinuity factors for flux reconstruction.

By default, the coupled x-y terms are neglected, as this requires knowledge of the average corner fluxes based on the adjacent nodes. Due to this restriction, the coefficients for the coupled x-y terms must be computed separately. Examine the source code in NEMDriver::fit_node_recon_params_corners for an example of how this should be done. It follows the ANOVA-HDMR Decomposition method [3].

Parameters:
  • dx (float) – Width of the node along the x axis.

  • dy (float) – Width of the node along the y axis.

  • keff (float) – Multiplication factor from reference Sn calculation.

  • xs (DiffusionCrossSection) – Diffusion cross sections homogenized for the node.

  • avg_flx (np.ndarray) – Average flux in each group within the node from reference Sn calculation.

  • j_x_neg (np.ndarray) – Reference net current in each group on the negative x boundary from the reference calculation.

  • j_x_pos (np.ndarray) – Reference net current in each group on the positive x boundary from the reference calculation.

  • j_y_neg (np.ndarray) – Reference net current in each group on the negative y boundary from the reference calculation.

  • j_y_pos (np.ndarray) – Reference net current in each group on the positive y boundary from the reference calculation.

dx#

Width along the x direction.

Type:

float

dy#

Width along the y direction.

Type:

float

keff#

Multiplication factor to use for the node.

Type:

float

ngroups#

Number of energy groups.

Type:

int

phi_0#

Average flux in the node for each group.

Type:

np.ndarray

flux_x#

Expansion of the nodal flux along the x axis.

Type:

NodalFlux1D

flux_y#

Expansion of the nodal flux along the y axis.

Type:

NodalFlux1D

__call__(x: float, y: float) ndarray[source]#

Evaluates the nodal flux in all groups at the position (x,y). The node is centered at (0,0), so the x value must be in the range [-dx/2, dx/2] and the y value must be in the range [-dy/2, dy/2].

Parameters:
  • x (float) – x position in the node.

  • y (float) – y position in the node.

Returns:

Nodal flux in all groups.

Return type:

np.ndarray

scarabee.reseau.compute_adf_cdf_from_cmfd(moc: MOCDriver, moc_to_nodal_cond_scheme: List[List[int]], symmetry: Symmetry, independent_quadrant: bool) Tuple[ndarray, ndarray][source]#

Computes the assembly and corner discontinuity factors in the few-group structure using the CMFD results. Currently, this method neglects the coupled x-y terms in the nodal construction when computing the CDFs. This simplifies the algorithm, but also ensures problems won’t occur if trying to run single-pin cell “assembly” calculations.

Parameters:
  • moc (MOCDriver) – Solved MOC problem for which the equivalence factors will be computed.

  • moc_to_nodal_cond_scheme (list of list of ints) – Scheme used to condense from the MOC problem to the resulting nodal diffusion parameters.

  • symmetry (Symmetry) – Indicates the symmetry of the problem.

  • independent_quadrant (bool) – If quarter symmetry is indicated (symmetry = Symmetry.Quarter), then unique ADFs and CDFs will still be generated for each side, as if full symmetry were being used.

Returns:

  • ADF (ndarray) – Assembly Discontinuity Factors

  • CDF (ndarray) – Corner Discontinuity Factors

Raises:

RuntimeError – If moc does not have a CMFD instance, or if it is not possible to create a condensation scheme to go from the CMFD group structure to the few-group structure.

scarabee.reseau.compute_adf_cdf_from_moc(moc: MOCDriver, moc_to_nodal_cond_scheme: List[List[int]], symmetry: Symmetry, independent_quadrant: bool) Tuple[ndarray, ndarray][source]#

Computes the assembly and corner discontinuity factors in the few-group structure using the MOC results.

Parameters:
  • moc (MOCDriver) – Solved MOC problem for which the equivalence factors will be computed.

  • moc_to_nodal_cond_scheme (list of list of ints) – Scheme used to condense from the MOC problem to the resulting nodal diffusion parameters.

  • symmetry (Symmetry) – Indicates the symmetry of the problem.

  • independent_quadrant (bool) – If quarter symmetry is indicated (symmetry = Symmetry.Quarter), then unique ADFs and CDFs will still be generated for each side, as if full symmetry were being used.

Returns:

  • ADF (ndarray) – Assembly Discontinuity Factors

  • CDF (ndarray) – Corner Discontinuity Factors

Full Core Tool Chains#

Building Full Core Models#

class scarabee.coeur.SimpleTile(diffusion_data: DiffusionData, form_factors: FormFactors)[source]#

Represents a single fuel assembly or reflector region, with diffusion cross sections, ADFs, CDFs, and form factors. Only a single set of cross sections are provided for the entire assembly. This is appropriate for symmetric fuel assemblies or reflector region.

Parameters:
  • diffusion_data (DiffusionData) – The diffusion cross sections, ADFs, and CDFs for the tile.

  • form_factors (FormFactors) – Form factors for the tile.

__init__(diffusion_data: DiffusionData, form_factors: FormFactors)[source]#
property num_x_slots: int#

Number of unique nodes which should be used along the x direction.

property num_y_slots: int#

Number of unique nodes which should be used along the y direction.

property x_width: float#

Width of the tile along the x direction.

property y_width: float#

Width of the tile along the y direction.

rotate_clockwise() CoreTile[source]#

Rotates the ADFs and form factors, corresponding to a 90 degree rotation of the assembly in the clockwise direction.

rotate_counterclockwise() CoreTile[source]#

Rotates the ADFs and form factors, corresponding to a 90 degree rotation of the assembly in the counter-clockwise direction.

reflect_across_x_axis() CoreTile[source]#

Performs a reflection of the ADFs and form factors across the x axis. This means that the +y and -y ADFs are swapped.

reflect_across_y_axis() CoreTile[source]#

Performs a reflection of the ADFs and form factors across the y axis. This means that the +x and -x ADFs are swapped.

class scarabee.coeur.QuadrantsTile(quad1: Optional[DiffusionData], quad2: Optional[DiffusionData], quad3: Optional[DiffusionData], quad4: Optional[DiffusionData], form_factors: FormFactors)[source]#

Represents a single fuel assembly with diffusion cross sections, ADFs, CDFs, and form factors. Each quadrant of the assembly has a unique DiffusionData instance. This is appropriate for asymmetric fuel assemblies. A single FormFactors instance is provided, which can be generated from independent quadrant assembly calculations.

Parameters:
  • quad1 (DiffusionData or None) – The diffusion cross sections, ADFs, and CDFs for the I quadrant.

  • quad2 (DiffusionData or None) – The diffusion cross sections, ADFs, and CDFs for the II quadrant.

  • quad3 (DiffusionData or None) – The diffusion cross sections, ADFs, and CDFs for the III quadrant.

  • quad4 (DiffusionData or None) – The diffusion cross sections, ADFs, and CDFs for the IV quadrant.

  • form_factors (FormFactors) – Form factors for the tile.

__init__(quad1: Optional[DiffusionData], quad2: Optional[DiffusionData], quad3: Optional[DiffusionData], quad4: Optional[DiffusionData], form_factors: FormFactors)[source]#
static from_independent_quadrant(diffusion_data: DiffusionData, form_factors: FormFactors) QuadrantsTile[source]#

Creates a QuadrantsTile from a single DiffusionData and FormFactors instance which were generated as independent quadrants, where ADFs and CDFs were generated along the symmetry axes of the assembly.

Parameters:
  • diffusion_data (DiffusionData) – Few-group cross sections, ADFs, and CDFs for the quarter assembly.

  • form_factors (FormFactors) – Form factors for the quarter assembly.

Return type:

QuadrantsTile

property num_x_slots: int#

Number of unique nodes which should be used along the x direction.

property num_y_slots: int#

Number of unique nodes which should be used along the y direction.

property x_width: float#

Width of the tile along the x direction.

property y_width: float#

Width of the tile along the y direction.

rotate_clockwise() CoreTile[source]#

Rotates the ADFs and form factors, corresponding to a 90 degree rotation of the assembly in the clockwise direction.

rotate_counterclockwise() CoreTile[source]#

Rotates the ADFs and form factors, corresponding to a 90 degree rotation of the assembly in the counter-clockwise direction.

reflect_across_x_axis() CoreTile[source]#

Performs a reflection of the ADFs and form factors across the x axis. This means that the +y and -y ADFs are swapped.

reflect_across_y_axis() CoreTile[source]#

Performs a reflection of the ADFs and form factors across the y axis. This means that the +x and -x ADFs are swapped.

class scarabee.coeur.CoreFormFactors(tile_width: float, num_tiles: int, form_factors: ndarray, z_widths: ndarray, symmetry: Symmetry = Symmetry.Full)[source]#

Evlauates the power form factors for any position across the entire core of the reactor. It is used to accomplish pin-power reconstruction after full core simulations.

Parameters:
  • tile_width (float) – Width of a single tile. Generally the assembly width/pitch.

  • num_tiles (int) – Number of tiles across the entire core at the widest point.

  • form_factors (3D array of FormFactors or 0.) – The array contains the form factors for all the assemblies in the core at all axial slices. It should be ordered in what you see is what you get (z,y,x), and the indexing for z and y is backwards (low z index means high z value). Entries for reflectors or zero flux regions should contain a 0.

  • z_widths (1D array of float) – Width of each axial slice from low to high z.

  • symmetry (Symmetry) – Symmetry used for the simulation. Default value is Symmetry.Full.

__init__(tile_width: float, num_tiles: int, form_factors: ndarray, z_widths: ndarray, symmetry: Symmetry = Symmetry.Full)[source]#
__call__(x, y, z)[source]#

Evaluates the form factors at the provided points.

Parameters:
  • x (float or np.ndarray of float) – x coordinate at which to evaluate the form factors

  • y (float or np.ndarray of float) – y coordinate at which to evaluate the form factors

  • z (float or np.ndarray of float) – z coordinate at which to evaluate the form factors

Returns:

If x, y, and z are all floats, then a single float with a form factor is returned. If any of them are an array, an array of form factors is returned.

Return type:

float or np.ndarray

class scarabee.coeur.CoreBuilder(tile_width: float, num_tiles: int, pitch: float, num_pins: int, tiles: ndarray, z_widths: ndarray, zmin_albedo: float, zmax_albedo: float)[source]#

The CoreBuilder class facilitates the construction of a NEMDiffusionDriver instance and a CoreFormFactors instnace, which are used to perform a full core simulation and compute assembly/pin powers. This class not not yet fully support defining cores with quarter symmetry in the x-y plane. Users should perform full core nodal simulations for the time being.

The user provided axial discretization is directly used in the nodal simulation, and will not be futher discretized.

Parameters:
  • tile_width (foat) – Width / pitch of an assembly in the core.

  • num_tiles (int) – Number of tiles along the diameter of the simulation. This is likely the number of assemblies across the core plus two (for a reflector tile on each side).

  • pitch (float) – The fuel pin pitch within an assembly.

  • num_pins (int) – Number of fuel pins along one side of an assembly.

  • tiles (np.ndarray of CoreTile or float) – A 3D array describing the geometry of the core. The array is indexed along z, y, then x, with z and y being indexed from large to small values (which permits a “what you see is what you get in your editor). Each entry in tiles should either be a SimpleTile or a QuadrantsTile. An entry of 0. can be used indicate a tile which is not used.

  • z_widths (np.ndarray) – A 2D array for the axial height of each axial slice of the geometry. This is the axial discretization which will be used in the nodal simulation.

  • zmin_albedo (float) – Albedo boundary condition at the -z boundary.

  • zmax_albedo (float) – Albedo boundary condition at the +z boundary.

__init__(tile_width: float, num_tiles: int, pitch: float, num_pins: int, tiles: ndarray, z_widths: ndarray, zmin_albedo: float, zmax_albedo: float)[source]#
compute_assembly_powers() ndarray[source]#

Computes the axially integrated average assembly powers. The powers are normalized such that the average power is unity, only considering tiles with a non-zero power.

Returns:

A 2D Numpy array with the average assembly powers.

Return type:

np.ndarray

compute_pin_powers(z: float | numpy.ndarray) ndarray[source]#

Computes the relative pin powers such that the average of all pin powers is unity, only considering pins with a non-zero power.

Parameters:

z (float or np.ndarray) – The z coordinates at which to evaluate the pin powers.

Returns:

A 3D Numpy array with the pin powers evaluated at the pin centers and the provided z points. Indexed according to x, y, and z.

Return type:

np.ndarray