openmc.UnstructuredMesh

class openmc.UnstructuredMesh(filename: str | PathLike, library: str, mesh_id: int | None = None, name: str = '', length_multiplier: float = 1.0, options: str | None = None)[source]

A 3D unstructured mesh

Added in version 0.12.

Changed in version 0.12.2: Support for libMesh unstructured meshes was added.

Parameters:
  • filename (path-like) – Location of the unstructured mesh file. Supported files for ‘moab’ library are .h5 and .vtk. Supported files for ‘libmesh’ library are exodus mesh files .exo.

  • library ({'moab', 'libmesh'}) – Mesh library used for the unstructured mesh tally

  • mesh_id (int) – Unique identifier for the mesh

  • name (str) – Name of the mesh

  • length_multiplier (float) – Constant multiplier to apply to mesh coordinates

  • options (str, optional) – Special options that control spatial search data structures used. This is currently only used to set parameters for MOAB’s AdaptiveKDTree. If None, OpenMC internally uses a default of “MAX_DEPTH=20;PLANE_SET=2;”.

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

  • name (str) – Name of the mesh

  • filename (str) – Name of the file containing the unstructured mesh

  • length_multiplier (float) – Multiplicative factor to apply to mesh coordinates

  • library ({'moab', 'libmesh'}) – Mesh library used for the unstructured mesh tally

  • options (str) –

    Special options that control spatial search data structures used. This is currently only used to set parameters for MOAB’s AdaptiveKDTree. If None, OpenMC internally uses a default of “MAX_DEPTH=20;PLANE_SET=2;”.

  • output (bool) – Indicates whether or not automatic tally output should be generated for this mesh

  • volumes (Iterable of float) – Volumes of the unstructured mesh elements

  • centroids (numpy.ndarray) – Centroids of the mesh elements with array shape (n_elements, 3)

  • vertices (numpy.ndarray) –

    Coordinates of the mesh vertices with array shape (n_elements, 3)

    Added in version 0.13.1.

  • connectivity (numpy.ndarray) –

    Connectivity of the elements with array shape (n_elements, 8)

    Added in version 0.13.1.

  • element_types (Iterable of integers) –

    Mesh element types

    Added in version 0.13.1.

  • total_volume (float) – Volume of the unstructured mesh in total

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

centroid(bin: int)[source]

Return the vertex averaged centroid of an element

Parameters:

bin (int) – Bin ID for the returned centroid

Returns:

x, y, z values of the element centroid

Return type:

numpy.ndarray

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 unstructured mesh object from XML element

Parameters:

elem (lxml.etree._Element) – XML element

Returns:

UnstructuredMesh generated from an XML element

Return type:

openmc.UnstructuredMesh

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

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.

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 volumes

Return Volumes for every mesh cell if populated by a StatePoint file

Returns:

volumes – Volumes

Return type:

numpy.ndarray

write_data_to_vtk(filename: str | PathLike | None = None, datasets: dict | None = None, volume_normalization: bool = True)[source]

Map data to unstructured VTK mesh elements.

If filename is None, then a filename will be generated based on the mesh ID, and exported to VTK format.

Parameters:
  • filename (str or pathlib.Path) – Name of the VTK file to write. If the filename ends in ‘.vtu’ then a binary VTU format file will be written, if the filename ends in ‘.vtk’ then a legacy VTK file will be written.

  • datasets (dict) – Dictionary whose keys are the data labels and values are numpy appropriately sized arrays of the data

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

write_vtk_mesh(**kwargs)[source]

Map data to unstructured VTK mesh elements.

Deprecated since version 0.13: Use UnstructuredMesh.write_data_to_vtk() instead.

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

  • datasets (dict) – Dictionary whose keys are the data labels and values are the data sets.

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