openmc.CylindricalMesh

class openmc.CylindricalMesh(r_grid: Sequence[float], z_grid: Sequence[float], phi_grid: Sequence[float] = (0, 6.283185307179586), origin: Sequence[float] = (0.0, 0.0, 0.0), mesh_id: int | None = None, name: str = '')[source]

A 3D cylindrical mesh

Parameters:
  • r_grid (numpy.ndarray) – 1-D array of mesh boundary points along the r-axis Requirement is r >= 0.

  • z_grid (numpy.ndarray) – 1-D array of mesh boundary points along the z-axis relative to the origin.

  • phi_grid (numpy.ndarray) – 1-D array of mesh boundary points along the phi-axis in radians. The default value is [0, 2π], i.e. the full phi range.

  • origin (numpy.ndarray) – 1-D array of length 3 the (x,y,z) origin of the mesh in cartesian coordinates

  • mesh_id (int) – Unique identifier for the mesh

  • name (str) – Name of the mesh

Variables:
  • id (int) – Unique identifier for the mesh

  • name (str) – Name of the mesh

  • dimension (Iterable of int) – The number of mesh cells in each direction (r_grid, phi_grid, z_grid).

  • n_dimension (int) – Number of mesh dimensions (always 3 for a CylindricalMesh).

  • r_grid (numpy.ndarray) – 1-D array of mesh boundary points along the r-axis. Requirement is r >= 0.

  • phi_grid (numpy.ndarray) – 1-D array of mesh boundary points along the phi-axis in radians. The default value is [0, 2π], i.e. the full phi range.

  • z_grid (numpy.ndarray) – 1-D array of mesh boundary points along the z-axis relative to the origin.

  • origin (numpy.ndarray) – 1-D array of length 3 the (x,y,z) origin of the mesh in cartesian coordinates

  • indices (Iterable of tuple) – An iterable of mesh indices for each mesh element, e.g. [(1, 1, 1), (2, 1, 1), …]

  • lower_left (Iterable of float) – The lower-left corner of the structured mesh. If only two coordinate are given, it is assumed that the mesh is an x-y mesh.

  • upper_right (Iterable of float) – The upper-right corner of the structured mesh. If only two coordinate are given, it is assumed that the mesh is an x-y mesh.

  • bounding_box (openmc.BoundingBox) – Axis-aligned bounding box of the mesh as defined by the upper-right and lower-left coordinates.

property centroids

Return coordinates of mesh element centroids.

Returns:

centroids – Returns a numpy.ndarray representing the mesh element centroid coordinates with a shape equal to (dim1, …, dimn, ndim). X, Y, Z values can be unpacked with xx, yy, zz = np.rollaxis(mesh.centroids, -1).

Return type:

numpy.ndarray

property centroids_cylindrical

Returns centroids of the mesh in cylindrical coordinates.

classmethod from_domain(domain: HasBoundingBox, dimension: Sequence[int] = (10, 10, 10), mesh_id: int | None = None, phi_grid_bounds: Sequence[float] = (0.0, 6.283185307179586), name: str = '', enclose_domain: bool = False)[source]

Create CylindricalMesh from a domain using its bounding box.

Parameters:
  • domain (HasBoundingBox) – The object passed in will be used as a template for this mesh. The bounding box of the property of the object passed will be used to set the r_grid, z_grid ranges.

  • dimension (Iterable of int) – The number of equally spaced mesh cells in each direction (r_grid, phi_grid, z_grid)

  • mesh_id (int) – Unique identifier for the mesh

  • phi_grid_bounds (numpy.ndarray) – Mesh bounds points along the phi-axis in radians. The default value is (0, 2π), i.e., the full phi range.

  • name (str) – Name of the mesh

  • enclose_domain (bool) – If True, the mesh will encompass the bounding box of the domain. If False, the mesh will be inscribed within the domain’s bounding box.

Returns:

CylindricalMesh instance

Return type:

openmc.CylindricalMesh

classmethod from_hdf5(group: Group, mesh_id: int, name: str)[source]

Create mesh from HDF5 group

Parameters:

group (h5py.Group) – Group in HDF5 file

Returns:

Instance of a MeshBase subclass

Return type:

openmc.MeshBase

classmethod from_xml_element(elem: Element)[source]

Generate a cylindrical mesh from an XML element

Parameters:

elem (lxml.etree._Element) – XML element

Returns:

Cylindrical mesh object

Return type:

openmc.CylindricalMesh

get_homogenized_materials(model: Model, n_samples: int | tuple[int, int, int] = 10000, include_void: bool = True, material_volumes: MeshMaterialVolumes | None = None, **kwargs) list[Material]

Generate homogenized materials over each element in a mesh.

Added in version 0.15.0.

Parameters:
  • model (openmc.Model) – Model containing materials to be homogenized and the associated geometry.

  • n_samples (int or 2-tuple of int) – Total number of rays to sample. The number of rays in each direction is determined by the aspect ratio of the mesh bounding box. When specified as a 3-tuple, it is interpreted as the number of rays in the x, y, and z dimensions.

  • include_void (bool, optional) – Whether homogenization should include voids.

  • material_volumes (MeshMaterialVolumes, optional) – Previously computed mesh material volumes to use for homogenization. If not provided, they will be computed by calling material_volumes().

  • **kwargs – Keyword-arguments passed to material_volumes().

Returns:

Homogenized material in each mesh element

Return type:

list of openmc.Material

get_indices_at_coords(coords: Sequence[float]) tuple[int, int, int][source]

Finds the index of the mesh voxel at the specified x,y,z coordinates.

Added in version 0.15.0.

Parameters:

coords (Sequence[float]) – The x, y, z axis coordinates

Returns:

The r, phi, z indices

Return type:

tuple[int, int, int]

material_volumes(model: Model, n_samples: int | tuple[int, int, int] = 10000, max_materials: int = 4, **kwargs) MeshMaterialVolumes

Determine volume of materials in each mesh element.

This method works by raytracing repeatedly through the mesh to count the estimated volume of each material in all mesh elements. Three sets of rays are used: one set parallel to the x-axis, one parallel to the y-axis, and one parallel to the z-axis.

Added in version 0.15.1.

Parameters:
  • model (openmc.Model) – Model containing materials.

  • n_samples (int or 3-tuple of int) – Total number of rays to sample. The number of rays in each direction is determined by the aspect ratio of the mesh bounding box. When specified as a 3-tuple, it is interpreted as the number of rays in the x, y, and z dimensions.

  • max_materials (int, optional) – Estimated maximum number of materials in any given mesh element.

  • **kwargs (dict) – Keyword arguments passed to openmc.lib.init()

Returns:

  • Dictionary-like object that maps material IDs to an array of volumes

  • equal in size to the number of mesh elements.

property midpoint_vertices

Create vertices that lie on the midpoint of element edges

classmethod reset_ids()

Reset counters

to_xml_element()[source]

Return XML representation of the mesh

Returns:

element – XML element containing mesh data

Return type:

lxml.etree._Element

property vertices
Return coordinates of mesh vertices in Cartesian coordinates. Also

see CylindricalMesh.vertices_cylindrical() and SphericalMesh.vertices_spherical() for coordinates in other coordinate systems.

Returns:

vertices – Returns a numpy.ndarray representing the coordinates of the mesh vertices with a shape equal to (dim1 + 1, …, dimn + 1, ndim). X, Y, Z values can be unpacked with xx, yy, zz = np.rollaxis(mesh.vertices, -1).

Return type:

numpy.ndarray

property vertices_cylindrical

Returns vertices of the mesh in cylindrical coordinates.

property volumes

Return Volumes for every mesh cell

Returns:

volumes – Volumes

Return type:

Iterable of float

write_data_to_vtk(filename: str | PathLike, datasets: dict | None = None, volume_normalization: bool = True, curvilinear: bool = False)

Creates a VTK object of the mesh

Parameters:
  • filename (str) – Name of the VTK file to write.

  • datasets (dict) – Dictionary whose keys are the data labels and values are the data sets. 1D datasets are expected to be extracted directly from statepoint data without reordering/reshaping. Multidimensional datasets are expected to have the same dimensions as the mesh itself with structured indexing in “C” ordering. See the “expand_dims” flag of get_reshaped_data() on reshaping tally data when using MeshFilter’s.

  • volume_normalization (bool, optional) – Whether or not to normalize the data by the volume of the mesh elements.

  • curvilinear (bool) – Whether or not to write curvilinear elements. Only applies to SphericalMesh and CylindricalMesh.

Raises:

ValueError – When the size of a dataset doesn’t match the number of mesh cells

Returns:

a VTK grid object representing the mesh

Return type:

vtk.StructuredGrid or vtk.UnstructuredGrid

Examples

1D data from a tally with only a mesh filter and heating score:

# pass the tally mean property of shape (N, 1, 1) directly to this # method; dimensions of size 1 will automatically removed >>> heating = tally.mean >>> mesh.write_data_to_vtk({‘heating’: heating})

Multidimensional data from a tally with only a mesh

# retrieve a data array with the mesh filter expanded into three # dimensions, ijk; additional dimensions of size one will # automatically be removed >>> heating = tally.get_reshaped_data(expand_dims=True) >>> mesh.write_data_to_vtk({‘heating’: heating})