armi.reactor.blocks module
Defines blocks, which are axial chunks of assemblies. They contain most of the state variables, including power, flux, and homogenized number densities.
Assemblies are made of blocks. Blocks are made of components.
- class armi.reactor.blocks.Block(name: str, height: float = 1.0)[source]
Bases:
Composite
A homogenized axial slab of material.
Blocks are stacked together to form assemblies.
Builds a new ARMI block.
- namestr
The name of this block
- heightfloat, optional
The height of the block in cm. Defaults to 1.0 so that
getVolume
assumes unit height.
- uniqID = 0
- PITCH_DIMENSION = 'op'
- pDefs = <armi.reactor.parameters.parameterDefinitions.ParameterDefinitionCollection object>
- createHomogenizedCopy(pinSpatialLocators=False)[source]
Create a copy of a block.
Notes
Used to implement a copy function for specific block types that can be much faster than a deepcopy by glossing over details that may be unnecessary in certain contexts.
This base class implementation is just a deepcopy of the block, in full detail (not homogenized).
- property core
- makeName(assemNum, axialIndex)[source]
Generate a standard block from assembly number.
This also sets the block-level assembly-num param.
Once, we used a axial-character suffix to represent the axial index, but this is inherently limited so we switched to a numerical name. The axial suffix needs can be brought in to plugins that require them.
Examples
>>> makeName(120, 5) 'B0120-005'
- getSmearDensity(cold=True)[source]
Compute the smear density of pins in this block.
Smear density is the area of the fuel divided by the area of the space available for fuel inside the cladding. Other space filled with solid materials is not considered available. If all the area is fuel, it has 100% smear density. Lower smear density allows more room for swelling.
Warning
This requires circular fuel and circular cladding. Designs that vary from this will be wrong. It may make sense in the future to put this somewhere a bit more design specific.
Notes
This only considers circular objects. If you have a cladding that is not a circle, it will be ignored.
Negative areas can exist for void gaps in the fuel pin. A negative area in a gap represents overlap area between two solid components. To account for this additional space within the pin cladding the abs(negativeArea) is added to the inner cladding area.
- autoCreateSpatialGrids(systemSpatialGrid=None)[source]
Creates a spatialGrid for a Block.
Blocks do not always have a spatialGrid from Blueprints, but some Blocks can have their spatialGrids inferred based on the multiplicty of their components. This would add the ability to create a spatialGrid for a Block and give its children the corresponding spatialLocators if certain conditions are met.
- Parameters:
systemSpatialGrid (Grid, optional) – Spatial Grid of the system-level parent of this Assembly that contains this Block.
- Raises:
ValueError – If the multiplicities of the block are not only 1 or N or if generated ringNumber leads to more positions than necessary.
- getMgFlux(adjoint=False, average=False, volume=None, gamma=False)[source]
Returns the multigroup neutron flux in [n/cm^2/s].
The first entry is the first energy group (fastest neutrons). Each additional group is the next energy group, as set in the ISOTXS library.
It is stored integrated over volume on self.p.mgFlux
- Parameters:
adjoint (bool, optional) – Return adjoint flux instead of real
average (bool, optional) – If true, will return average flux between latest and previous. Doesn’t work for pin detailed yet
volume (float, optional) – If average=True, the volume-integrated flux is divided by volume before being returned. The user may specify a volume, or the function will obtain the block volume directly.
gamma (bool, optional) – Whether to return the neutron flux or the gamma flux.
- Returns:
flux
- Return type:
multigroup neutron flux in [n/cm^2/s]
- setPinMgFluxes(fluxes, adjoint=False, gamma=False)[source]
Store the pin-detailed multi-group neutron flux.
- Parameters:
fluxes (np.ndarray) – The block-level pin multigroup fluxes. fluxes[i, g] represents the flux in group g for pin i. Flux units are the standard n/cm^2/s. The “ARMI pin ordering” is used, which is counter-clockwise from 3 o’clock.
adjoint (bool, optional) – Whether to set real or adjoint data.
gamma (bool, optional) – Whether to set gamma or neutron data.
Outputs –
------- –
self.p.pinMgFluxes (np.ndarray) – The block-level pin multigroup fluxes. pinMgFluxes[i, g] represents the flux in group g for pin i. Flux units are the standard n/cm^2/s. The “ARMI pin ordering” is used, which is counter-clockwise from 3 o’clock.
- getMicroSuffix()[source]
Returns the microscopic library suffix (e.g. ‘AB’) for this block.
DIF3D and MC2 are limited to 6 character nuclide labels. ARMI by convention uses the first 4 for nuclide name (e.g. U235, PU39, etc.) and then uses the 5th character for cross-section type and the 6th for burnup group. This allows a variety of XS sets to be built modeling substantially different blocks.
Notes
The single-letter use for xsType and envGroup limit users to 52 groups of each. ARMI will allow 2-letter xsType designations if and only if the envGroup setting has length 1 (i.e. no burnup/temp groups are defined). This is useful for high-fidelity XS modeling of V&V models such as the ZPPRs.
- setHeight(modifiedHeight, conserveMass=False, adjustList=None)[source]
Set a new height of the block.
- Parameters:
Notes
There is a coupling between block heights, the parent assembly axial mesh, and the ztop/zbottom/z params of the sibling blocks. When you set a height, all those things are invalidated. Thus, this method has to go through and update them via
parent.calculateZCoords
. This could be inefficient though it has not been identified as a bottleneck. Possible improvements include deriving z/ztop/zbottom on the fly and invalidating the parent mesh with some kind of flag, signaling it to recompute itself on demand. Developers can get around some of the O(N^2) scaling of this by settingp.height
directly but they must know to update the dependent objects after they do that. Use with care.See also
armi.reactor.reactors.Core.updateAxialMesh
May need to be called after this.
armi.reactor.assemblies.Assembly.calculateZCoords
Recalculates z-coords, automatically called by this.
- getFlowAreaPerPin()[source]
Return the flowing coolant area of the block in cm^2, normalized to the number of pins in the block.
NumPins looks for max number of fuel, clad, control, etc.
See also
armi.reactor.blocks.Block.getNumPins
figures out numPins
- adjustUEnrich(newEnrich)[source]
Adjust U-235/U-238 mass ratio to a mass enrichment.
- Parameters:
newEnrich (float) – New U-235 enrichment in mass fraction
Notes
completeInitialLoading must be run because adjusting the enrichment actually changes the mass slightly and you can get negative burnups, which you do not want.
- getArea(cold=False)[source]
Return the area of a block for a full core or a 1/3 core model.
Area is consistent with the area in the model, so if you have a central assembly in a 1/3 symmetric model, this will return 1/3 of the total area of the physical assembly. This way, if you take the sum of the areas in the core (or count the atoms in the core, etc.), you will have the proper number after multiplying by the model symmetry.
- Parameters:
cold (bool) – flag to indicate that cold (as input) dimensions are required
Notes
This might not work for a 1/6 core model (due to symmetry line issues).
- Returns:
area
- Return type:
float (cm^2)
See also
armi.reactor.blocks.Block.getMaxArea
return the full area of the physical assembly disregarding model symmetry
- getVolume()[source]
Return the volume of a block.
- Returns:
volume – Block or component volume in cm^3
- Return type:
- getSymmetryFactor()[source]
Return a scaling factor due to symmetry on the area of the block or its components.
Takes into account assemblies that are bisected or trisected by symmetry lines
In 1/3 symmetric cases, the central assembly is 1/3 a full area. If edge assemblies are included in a model, the symmetry factor along both edges for overhanging assemblies should be 2.0. However, ARMI runs in most scenarios with those assemblies on the 120-edge removed, so the symmetry factor should generally be just 1.0.
See also
armi.reactor.converters.geometryConverter.EdgeAssemblyChanger.scaleParamsRelatedToSymmetry
- adjustDensity(frac, adjustList, returnMass=False)[source]
adjusts the total density of each nuclide in adjustList by frac.
- Parameters:
- Returns:
mass – Mass difference in grams. If you subtract mass, mass will be negative. If returnMass is False (default), this will always be zero.
- Return type:
- completeInitialLoading(bolBlock=None)[source]
Does some BOL bookkeeping to track things like BOL HM density for burnup tracking.
This should run after this block is loaded up at BOC (called from Reactor.initialLoading).
The original purpose of this was to get the moles HM at BOC for the moles Pu/moles HM at BOL calculation.
This also must be called after modifying something like the smear density or zr fraction in an optimization case. In ECPT cases, a BOL block must be passed or else the burnup will try to get based on a pre-burned value.
- Parameters:
bolBlock (Block, optional) – A BOL-state block of this block type, required for perturbed equilibrium cases. Must have the same enrichment as this block!
- Returns:
hmDens – The heavy metal number density of this block.
- Return type:
See also
Reactor.importGeom
,depletion._updateBlockParametersAfterDepletion
- setB10VolParam(heightHot)[source]
Set the b.p.initialB10ComponentVol param according to the volume of boron-10 containing components.
- Parameters:
heightHot (Boolean) – True if self.height() is cold height
- replaceBlockWithBlock(bReplacement)[source]
Replace the current block with the replacementBlock.
Typically used in the insertion of control rods.
- getComponentsThatAreLinkedTo(comp, dim)[source]
Determine which dimensions of which components are linked to a specific dimension of a particular component.
Useful for breaking fuel components up into individuals and making sure anything that was linked to the fuel mult (like the cladding mult) stays correct.
- getComponentsInLinkedOrder(componentList=None)[source]
Return a list of the components in order of their linked-dimension dependencies.
- Parameters:
components (list, optional) – A list of components to consider. If None, this block’s components will be used.
Notes
This means that components other components are linked to come first.
- getSortedComponentsInsideOfComponent(component)[source]
Returns a list of components inside of the given component sorted from innermost to outermost.
- Parameters:
component (object) – Component to look inside of.
Notes
If you just want sorted components in this block, use
sorted(self)
. This will never include anyDerivedShape
objects. Since they have a derived area they don’t have a well-defined dimension. For now we just ignore them. If they are desired in the future some knowledge of their dimension will be required while they are being derived.
- mergeWithBlock(otherBlock, fraction)[source]
Turns this block into a mixture of this block and some other block.
- Parameters:
Notes
This merges on a high level (using number densities). Components will not be merged.
This is used e.g. for inserting a control block partially to get a very tight criticality control. In this case, a control block would be merged with a duct block. It is also used when a control rod is specified as a certain length but that length does not fit exactly into a full block.
- getComponentAreaFrac(typeSpec)[source]
Returns the area fraction of the specified component(s) among all components in the block.
- Parameters:
typeSpec (Flags or list of Flags) – Component types to look up
Examples
>>> b.getComponentAreaFrac(Flags.CLAD) 0.15
- Returns:
The area fraction of the component.
- Return type:
- getDim(typeSpec, dimName)[source]
Search through blocks in this assembly and find the first component of compName. Then, look on that component for dimName.
- Parameters:
- Returns:
dimVal – The dimension in cm.
- Return type:
Examples
>>> getDim(Flags.WIRE,'od') 0.01
- getPinCenterFlatToFlat(cold=False)[source]
Return the flat-to-flat distance between the centers of opposing pins in the outermost ring.
- getPitch(returnComp=False)[source]
Return the center-to-center hex pitch of this block.
- Parameters:
returnComp (bool, optional) – If true, will return the component that has the maximum pitch as well
- Returns:
pitch (float or None) – Hex pitch in cm, if well-defined. If there is no clear component for determining pitch, returns None
component (Component or None) – Component that has the max pitch, if returnComp == True. If no component is found to define the pitch, returns None
Notes
The block stores a reference to the component that defines the pitch, making the assumption that while the dimensions can change, the component containing the largest dimension will not. This lets us skip the search for largest component. We still need to ask the largest component for its current dimension in case its temperature changed, or was otherwise modified.
See also
setPitch
sets pitch
- getPinPitch(cold=False)[source]
Return sub-block pitch in blocks.
This assumes the spatial grid is defined by unit steps
- getLargestComponent(dimension)[source]
Find the component with the largest dimension of the specified type.
- Parameters:
dimension (str) – The name of the dimension to find the largest component of.
- Returns:
largestComponent – The component with the largest dimension of the specified type.
- Return type:
armi.reactor.components.Component
- setPitch(val, updateBolParams=False)[source]
Sets outer pitch to some new value.
This sets the settingPitch and actually sets the dimension of the outer hexagon.
During a load (importGeom), the setDimension doesn’t usually do anything except set the setting See Issue 034
But during a actual case modification (e.g. in an optimization sweep, then the dimension has to be set as well.
See also
getPitch
gets the pitch
- getMfp(gamma=False)[source]
Calculate the mean free path for neutron or gammas in this block.
\[<\Sigma> = \frac{\sum_E(\phi_e \Sigma_e dE)}{\sum_E (\phi_e dE)} = \frac{\sum_E(\phi_e N \sum_{\text{type}}(\sigma_e) dE}{\sum_E (\phi_e dE))}\]Block macro is the sum of macros of all nuclides.
phi_g = flux*dE already in multigroup method.
- getBlocks()[source]
This method returns all the block(s) included in this block its implemented so that methods could iterate over reactors, assemblies or single blocks without checking to see what the type of the reactor-family object is.
- updateComponentDims()[source]
This method updates all the dimensions of the components.
Notes
This is VERY useful for defining a ThRZ core out of differentialRadialSegements whose dimensions are connected together some of these dimensions are derivative and can be updated by changing dimensions in a Parameter Component or other linked components
See also
armi.reactor.components.DifferentialRadialSegment.updateDims
,armi.reactor.components.Parameters
,armi.physics.optimize.OptimizationInterface.modifyCase
- breakFuelComponentsIntoIndividuals()[source]
Split block-level components (in fuel blocks) into pin-level components.
The fuel component will be broken up according to its multiplicity.
Order matters! The first pin component will be located at a particular (x, y), which will be used in the fluxRecon module to determine the interpolated flux.
The fuel will become fuel001 through fuel169 if there are 169 pins.
- getIntegratedMgFlux(adjoint=False, gamma=False)[source]
Return the volume integrated multigroup neutron tracklength in [n-cm/s].
The first entry is the first energy group (fastest neutrons). Each additional group is the next energy group, as set in the ISOTXS library.
- getLumpedFissionProductCollection()[source]
Get collection of LFP objects. Will work for global or block-level LFP models.
- Returns:
lfps – lfpName keys , lfp object values
- Return type:
LumpedFissionProduct
- rotate(rad)[source]
Function for rotating a block’s spatially varying variables by a specified angle (radians).
- Parameters:
rad (float) – Number (in radians) specifying the angle of counter clockwise rotation.
- setAxialExpTargetComp(targetComponent)[source]
Sets the targetComponent for the axial expansion changer.
Parameter
- targetComponent:
Component
object Component specified to be target component for axial expansion changer
- targetComponent:
- getPinLocations() list[armi.reactor.grids.locations.IndexLocation] [source]
Produce all the index locations for pins in the block.
- Returns:
Integer locations where pins can be found in the block.
- Return type:
list[grids.IndexLocation]
Notes
Only components with
Flags.CLAD
are considered to define a pin’s location.See also
- getPinCoordinates() ndarray [source]
Compute the local centroid coordinates of any pins in this block.
The pins must have a CLAD-flagged component for this to work.
- Returns:
localCoords –
(N, 3)
array of coordinates for pins locations.localCoords[i]
contains a triplet of the x, y, z location for pini
. Ordered according to how they are listed as children- Return type:
numpy.ndarray
See also
- getTotalEnergyGenerationConstants()[source]
Get the total energy generation group constants for a block.
Gives the total energy generation rates when multiplied by the multigroup flux.
- Returns:
totalEnergyGenConstant – Total (fission + capture) energy generation group constants (Joules/cm)
- Return type:
np.ndarray
- getFissionEnergyGenerationConstants()[source]
Get the fission energy generation group constants for a block.
Gives the fission energy generation rates when multiplied by the multigroup flux.
- Returns:
fissionEnergyGenConstant – Energy generation group constants (Joules/cm)
- Return type:
np.ndarray
- Raises:
RuntimeError: – Reports if a cross section library is not assigned to a reactor.
- getCaptureEnergyGenerationConstants()[source]
Get the capture energy generation group constants for a block.
Gives the capture energy generation rates when multiplied by the multigroup flux.
- Returns:
fissionEnergyGenConstant – Energy generation group constants (Joules/cm)
- Return type:
np.ndarray
- Raises:
RuntimeError: – Reports if a cross section library is not assigned to a reactor.
- getNeutronEnergyDepositionConstants()[source]
Get the neutron energy deposition group constants for a block.
- Returns:
energyDepConstants – Neutron energy generation group constants (in Joules/cm)
- Return type:
np.ndarray
- Raises:
RuntimeError: – Reports if a cross section library is not assigned to a reactor.
- getGammaEnergyDepositionConstants()[source]
Get the gamma energy deposition group constants for a block.
- Returns:
energyDepConstants – Energy generation group constants (in Joules/cm)
- Return type:
np.ndarray
- Raises:
RuntimeError: – Reports if a cross section library is not assigned to a reactor.
- getInputHeight() float [source]
Determine the input height from blueprints.
- Returns:
Height for this block pulled from the blueprints.
- Return type:
- Raises:
AttributeError – If no ancestor of this block contains the input blueprints. Blueprints are usually stored on the reactor object, which is typically an ancestor of the block (block -> assembly -> core -> reactor). However, this may be the case when creating blocks from scratch in testing where the entire composite tree may not exist.
- paramCollectionType
alias of
BlockParameterCollection
- class armi.reactor.blocks.HexBlock(name, height=1.0)[source]
Bases:
Block
Defines a HexBlock.
- PITCH_COMPONENT_TYPE: ClassVar[Optional[Tuple[Type[Component], ...]]] = (<class 'armi.reactor.components.basicShapes.Hexagon'>,)
- createHomogenizedCopy(pinSpatialLocators=False)[source]
Create a new homogenized copy of a block that is less expensive than a full deepcopy.
Notes
This can be used to improve performance when a new copy of a reactor needs to be built, but the full detail of the block (including component geometry, material, number density, etc.) is not required for the targeted physics solver being applied to the new reactor model.
The main use case is for the uniform mesh converter (UMC). Frequently, a deterministic neutronics solver will require a uniform mesh reactor, which is produced by the UMC. Many deterministic solvers for fast spectrum reactors will also treat the individual blocks as homogenized mixtures. Since the neutronics solver does not need to know about the geometric and material details of the individual child components within a block, we can save significant effort while building the uniform mesh reactor with the UMC by omitting this detailed data and only providing the necessary level of detail for the uniform mesh reactor: number densities on each block.
Individual components within a block can have different temperatures, and this can affect cross sections. This temperature variation is captured by the lattice physics module. As long as temperature distribution is correctly captured during cross section generation, it doesn’t need to be transferred to the neutronics solver directly through this copy operation.
If you make a new block, you must add it to an assembly and a reactor.
- Returns:
A homogenized block containing a single Hexagon Component that contains an average temperature and the number densities from the original block.
- Return type:
b
- setPinPowers(powers, powerKeySuffix='')[source]
Updates the pin linear power densities of this block for the current rotation. The linear densities are represented by the linPowByPin parameter.
It is assumed that
initializePinLocations()
has already been executed for fueled blocks in order to access the pinLocation parameter. The pinLocation parameter is not accessed for non-fueled blocks.The linPowByPin parameter can be directly assigned to instead of using this method if the multiplicity of the pins in the block is equal to the number of pins in the block.
- Parameters:
powers (list of floats, required) – The block-level pin linear power densities. powers[i] represents the average linear power density of pin i. The units of linear power density is watts/cm (i.e., watts produced per cm of pin length). The “ARMI pin ordering” must be be used, which is counter-clockwise from 3 o’clock.
powerKeySuffix (str, optional) – Must be either an empty string,
NEUTRON
, orGAMMA
. Defaults to empty string.
Notes
This method can handle assembly rotations by using the pinLocation parameter.
- rotate(rad: float)[source]
Rotates a block’s spatially varying parameters by a specified angle in the counter-clockwise direction.
The parameters must have a ParamLocation of either CORNERS or EDGES and must be a Python list of length 6 in order to be eligible for rotation; all parameters that do not meet these two criteria are not rotated.
- Parameters:
rad (float, required) – Angle of counter-clockwise rotation in units of radians. Rotations must be in 60-degree increments (i.e., PI/6, PI/3, PI, 2 * PI/3, 5 * PI/6, and 2 * PI)
- getPinToDuctGap(cold=False)[source]
Returns the distance in cm between the outer most pin and the duct in a block.
- Parameters:
cold (boolean) – Determines whether the results should be cold or hot dimensions.
- Returns:
pinToDuctGap – Returns the diameteral gap between the outer most pins in a hex pack to the duct inner face to face in cm.
- Return type:
- getRotationNum() int [source]
Get index 0 through 5 indicating number of rotations counterclockwise around the z-axis.
- setRotationNum(rotNum: int)[source]
Set orientation based on a number 0 through 5 indicating number of rotations counterclockwise around the z-axis.
- getSymmetryFactor()[source]
Return a factor between 1 and N where 1/N is how much cut-off by symmetry lines this mesh cell is.
Reactor-level meshes have symmetry information so we have a reactor for this to work. That’s why it’s not implemented on the grid/locator level.
When edge-assemblies are included on both edges (i.e. MCNP or DIF3D-FD 1/3-symmetric cases), the edge assemblies have symmetry factors of 2.0. Otherwise (DIF3D-nodal) there’s a full assembly on the bottom edge (overhanging) and no assembly at the top edge so the ones at the bottom are considered full (symmetryFactor=1).
If this block is not in any grid at all, then there can be no symmetry so return 1.
- autoCreateSpatialGrids(systemSpatialGrid=None)[source]
Given a block without a spatialGrid, create a spatialGrid and give its children the corresponding spatialLocators (if it is a simple block).
In this case, a simple block would be one that has either multiplicity of components equal to 1 or N but no other multiplicities. Also, this should only happen when N fits exactly into a given number of hex rings. Otherwise, do not create a grid for this block.
- Parameters:
systemSpatialGrid (Grid, optional) – Spatial Grid of the system-level parent of this Assembly that contains this Block.
Notes
When a hex grid has another hex grid nested inside it, the nested grid has the opposite orientation (corners vs flats up). This method takes care of that.
If components inside this block are multiplicity 1, they get a single locator at the center of the grid cell. If the multiplicity is greater than 1, all the components are added to a multiIndexLocation on the hex grid.
- Raises:
ValueError – If the multiplicities of the block are not only 1 or N or if generated ringNumber leads to more positions than necessary.
- getPinCenterFlatToFlat(cold=False)[source]
Return the flat-to-flat distance between the centers of opposing pins in the outermost ring.
- getPinPitch(cold=False)[source]
Get the pin pitch in cm.
Assumes that the pin pitch is defined entirely by contacting cladding tubes and wire wraps. Grid spacers not yet supported.
- Parameters:
cold (boolean) – Determines whether the dimensions should be cold or hot
- Returns:
pinPitch – pin pitch in cm
- Return type:
- getHydraulicDiameter()[source]
Return the hydraulic diameter in this block in cm.
Hydraulic diameter is 4A/P where A is the flow area and P is the wetted perimeter. In a hex assembly, the wetted perimeter includes the cladding, the wire wrap, and the inside of the duct. The flow area is the inner area of the duct minus the area of the pins and the wire.
- paramCollectionType
alias of
BlockParameterCollection
- class armi.reactor.blocks.CartesianBlock(name: str, height: float = 1.0)[source]
Bases:
Block
Builds a new ARMI block.
- namestr
The name of this block
- heightfloat, optional
The height of the block in cm. Defaults to 1.0 so that
getVolume
assumes unit height.
- PITCH_DIMENSION = 'widthOuter'
- getSymmetryFactor()[source]
Return a factor between 1 and N where 1/N is how much cut-off by symmetry lines this mesh cell is.
- getPinCenterFlatToFlat(cold=False)[source]
Return the flat-to-flat distance between the centers of opposing pins in the outermost ring.
- paramCollectionType
alias of
BlockParameterCollection
- class armi.reactor.blocks.ThRZBlock(name: str, height: float = 1.0)[source]
Bases:
Block
Builds a new ARMI block.
- namestr
The name of this block
- heightfloat, optional
The height of the block in cm. Defaults to 1.0 so that
getVolume
assumes unit height.
- paramCollectionType
alias of
BlockParameterCollection