openmc.RectLattice¶
- class openmc.RectLattice(lattice_id=None, name='')[source]¶
A lattice consisting of rectangular prisms.
To completely define a rectangular lattice, the
RectLattice.lower_left
RectLattice.pitch
,RectLattice.outer
, andRectLattice.universes
properties need to be set.Most methods for this class use a natural indexing scheme wherein elements are assigned an index corresponding to their position relative to the (x,y,z) axes in a Cartesian coordinate system, i.e., an index of (0,0,0) in the lattice gives the element whose x, y, and z coordinates are the smallest. However, note that when universes are assigned to lattice elements using the
RectLattice.universes
property, the array indices do not correspond to natural indices.- Parameters
- Variables
id (int) – Unique identifier for the lattice
name (str) – Name of the lattice
pitch (Iterable of float) – Pitch of the lattice in the x, y, and (if applicable) z directions in cm.
outer (openmc.Universe) – A universe to fill all space outside the lattice
universes (Iterable of Iterable of openmc.Universe) – A two- or three-dimensional list/array of universes filling each element of the lattice. The first dimension corresponds to the z-direction (if applicable), the second dimension corresponds to the y-direction, and the third dimension corresponds to the x-direction. Note that for the y-direction, a higher index corresponds to a lower physical y-value. Each z-slice in the array can be thought of as a top-down view of the lattice.
lower_left (Iterable of float) – The Cartesian coordinates of the lower-left corner of the lattice. If the lattice is two-dimensional, only the x- and y-coordinates are specified.
indices (list of tuple) – A list of all possible (z,y,x) or (y,x) lattice element indices. These indices correspond to indices in the
RectLattice.universes
property.ndim (int) – The number of dimensions of the lattice
shape (Iterable of int) – An array of two or three integers representing the number of lattice cells in the x- and y- (and z-) directions, respectively.
- create_xml_subelement(xml_element, memo=None)[source]¶
Add the lattice xml representation to an incoming xml element
- Parameters
xml_element (xml.etree.ElementTree.Element) – XML element to be added to
memo (set or None) – A set of object id’s representing geometry entities already written to the xml_element. This parameter is used internally and should not be specified by users.
- Returns
- Return type
- discretize(strategy='degenerate', universes_to_ignore=[], materials_to_clone=[], lattice_neighbors=[], key=<function RectLattice.<lambda>>)[source]¶
Discretize the lattice with either a degenerate or a local neighbor symmetry strategy
‘Degenerate’ clones every universe in the lattice, thus making them all uniquely defined. This is typically required if depletion or thermal hydraulics will make every universe’s environment unique.
‘Local neighbor symmetry’ groups universes with similar neighborhoods. These clusters of cells and materials provide increased convergence speed to multi-group cross sections tallies. The user can specify the lattice’s neighbors to discriminate between two sides of a lattice for example.
- Parameters
strategy ({'degenerate', 'lns'}) – Which strategy to adopt when discretizing the lattice
universes_to_ignore (Iterable of Universe) – Lattice universes that need not be discretized
materials_to_clone (Iterable of Material) – List of materials that should be cloned when discretizing
lattice_neighbors (Iterable of Universe) – List of the lattice’s neighbors. By default, if present, the lattice outer universe will be used. The neighbors should be ordered as follows [top left, top, top right, left, right, bottom left, bottom, bottom right]
key (function) – Function of argument a universe that is used to extract a comparison key. This function will be called on each universe’s neighbors in the lattice to form a neighbor pattern. This pattern is then used to identify unique neighbor symmetries.
- find_element(point)[source]¶
Determine index of lattice element and local coordinates for a point
- Parameters
point (Iterable of float) – Cartesian coordinates of point
- Returns
2- or 3-tuple of int – A tuple of the corresponding (x,y,z) lattice element indices
3-tuple of float – Carestian coordinates of the point in the corresponding lattice element coordinate system
- classmethod from_hdf5(group, universes)[source]¶
Create rectangular lattice from HDF5 group
- Parameters
group (h5py.Group) – Group in HDF5 file
universes (dict) – Dictionary mapping universe IDs to instances of
openmc.Universe
.
- Returns
Rectangular lattice
- Return type
- classmethod from_xml_element(elem, get_universe)[source]¶
Generate rectangular lattice from XML element
- Parameters
elem (xml.etree.ElementTree.Element) – <lattice> element
get_universe (function) – Function returning universe (defined in
openmc.Geometry.from_xml()
)
- Returns
Rectangular lattice
- Return type
- get_local_coordinates(point, idx)[source]¶
Determine local coordinates of a point within a lattice element
- Parameters
point (Iterable of float) – Cartesian coordinates of point
idx (Iterable of int) – (x,y,z) indices of lattice element. If the lattice is 2D, the z index can be omitted.
- Returns
Cartesian coordinates of point in the lattice element coordinate system
- Return type
3-tuple of float
- get_universe_index(idx)[source]¶
Return index in the universes array corresponding to a lattice element index
- Parameters
idx (Iterable of int) – Lattice element indices in the \((x,y,z)\) coordinate system
- Returns
Indices used when setting the
RectLattice.universes
property- Return type
2- or 3-tuple of int