Library design & HOWTO

Memory and file management

The library is designed to give you control over memory and file management: you have to explicitly open the data files and pass these opened files to the library to get back NumPy arrays with the data. This adds some bookkeeping overhead for you, but it enables you to choose between:

  1. traditional I/O and getting arrays which own a copy of the data, or

  2. memory-mapped I/O and getting arrays that are views into the binary EnSight files on disk.

Depending on your use-case, using memory-mapped I/O can be very efficient. It also allows for limited editing of EnSight files (eg. transforming coordinates or variable values in-place).

The following code example demonstrates this feature:

import ensightreader

case = ensightreader.read_case("data/sphere/sphere.case")
geofile = case.get_geometry_model()
part_ids = geofile.get_part_ids()
part = geofile.get_part_by_id(part_ids[0])

# (1) traditional I/O
with geofile.open() as fp_geo:
    nodes = part.read_nodes(fp_geo)  # owned buffer

# (2.1) read-only memory-mapped I/O
with geofile.mmap() as mm_geo:
    nodes = part.read_nodes(mm_geo)  # buffer backed by read-only mmap

# (2.2) write-through memory-mapped I/O
with geofile.mmap_writable() as mm_geo:
    nodes = part.read_nodes(mm_geo)  # buffer backed by write-through mmap
    nodes[:, 0] = 0.0                # set X coordinate to zero for part nodes

See also

The ensight_transform.py and ensight2vtk.py scripts for examples how to use memory-mapped I/O with the library. It’s done by simply passing mmap objects wherever the library expects opened files. Make sure to keep the mmap objects opened for as long as you need the arrays that are read from them!

Comparison with VTK library

The Visualisation Toolkit (VTK) is a great library for handling scientific and engineering data. Among its many features it comes with a reader (and a somewhat less featureful writer) for the EnSight Gold format.

Both libraries are quite different in focus. Put shortly, ensight-reader is a library you would use if you wanted to implement an EnSight Gold reader in a high-level library such as VTK.

Common points

  • both can be used from Python (VTK is a C++ library with official Python bindings)

  • both give you access to node coordinates, cell connectivity, variable data as Numpy arrays

ensight-reader specifics

  • supports partial reading (you can read just one variable for one part, etc.)

  • supports memory-mapped access (to minimize copying)

  • doesn’t do any conversion/interpretation of the data (useful if you’re doing something low-level yourself)

VTK specifics

  • supports some features that ensight-reader currently does not (eg. structured parts)

  • has more high-level API, can do much more with the data (eg. interpolate node/element data, save to different formats, show in 3D viewport, …)

  • does not support partial reading of parts

  • does not support memory-mapped access

  • converts the data into VTK datastructures (vtkUnstructuredGrid, etc.)

Code example

VTK library

>>> import vtk
>>> from vtk.numpy_interface import dataset_adapter as dsa

>>> reader = vtk.vtkEnSightGoldBinaryReader()
>>> reader.SetCaseFileName("data/sphere/sphere.case")
>>> reader.Update()
>>> case = reader.GetOutput()

>>> part = case.GetBlock(0)
>>> part_ = dsa.WrapDataObject(part)
>>> part_.Points
VTKArray([[ 0.0000000e+00,  0.0000000e+00,  5.0000000e+00],
          [ 0.0000000e+00,  0.0000000e+00, -5.0000000e+00],
          ...
          [ 1.5340106e+00, -1.5340106e+00, -4.5048442e+00]], dtype=float32)

>>> part_.PointData["RTData"]
VTKArray([220.84135, 220.84135, 223.80856, 233.50835, 217.5993 ,
          ...
          213.36838, 210.3635 , 210.3635 , 213.36838, 232.34589],
         dtype=float32)

>>> part_.Cells
VTKArray([ 3,  2,  8,  0,  3,  8, 14,  0,  3, 14, 20,  0,  3, 20, 26,  0,
          ...
           3, 47, 48,  6,  3, 47,  6,  5,  3, 48, 49,  7,  3, 48,  7,  6],
         dtype=int64)

ensight-reader library

>>> import ensightreader

>>> case = ensightreader.read_case("data/sphere/sphere.case")
>>> geofile = case.get_geometry_model()
>>> part_ids = geofile.get_part_ids()
>>> part = geofile.get_part_by_id(part_ids[0])

>>> with geofile.open() as fp_geo:
...     nodes = part.read_nodes(fp_geo)
>>> nodes
array([[ 0.0000000e+00,  0.0000000e+00,  5.0000000e+00],
       ...
       [ 1.5340106e+00, -1.5340106e+00, -4.5048442e+00]], dtype=float32)

>>> with geofile.open() as fp_geo:
...     block = part.element_blocks[0]
...     connectivity = block.read_connectivity(fp_geo)
>>> connectivity
array([[ 3,  9,  1],
       ...
       [49,  8,  7]])

>>> variable = case.get_variable("RTData")
>>> with variable.open() as fp_var:
...     variable_data = variable.read_node_data(fp_var, part.part_id)
>>> variable_data
array([220.84135, 220.84135, 223.80856, 233.50835, 217.5993 , 217.5993 ,
       ...
       213.36838, 232.34589], dtype=float32)

Writing data using ensight-reader

Despite its name, ensight-reader can be used to modify existing EnSight Gold cases or even create cases from scratch. Due to the unopinionated, low-level nature of the library, this can be tricky – especially if you’re not already familiar with inner workings of the EnSight Gold format.

Here are some examples to get you started:

Modifying variable values

>>> import ensightreader

>>> case = ensightreader.read_case("data/sphere/sphere.case")
>>> geofile = case.get_geometry_model()
>>> part_ids = geofile.get_part_ids()

>>> if case.variables["RTData"].timeset is None:
...     timesteps = [0]
... else:
...     timesteps = list(range(len(case.get_time_values()))

>>> for timestep in timesteps:
...     variable = case.get_variable("RTData", timestep)
...     with variable.mmap_writable() as mm_var:
...         for part_id in part_ids:
...             arr = variable.read_node_data(fp_var, part_id)
...             arr[:] *= 2

Defining new per-case constant

>>> import ensightreader

>>> case = ensightreader.read_case("data/sphere/sphere.case")
>>> constant = ensightreader.EnsightConstantVariable(timeset=None, variable_name="my_variable", values=[42.0])
>>> case.constant_variables[constant.variable_name] = constant
>>> case.to_file("data/sphere/sphere-with-constant.case")

Changing node coordinates

>>> import ensightreader

>>> case = ensightreader.read_case("data/sphere/sphere.case")
>>> geofile = case.get_geometry_model()
>>> part_ids = geofile.get_part_ids()

>>> with geofile.mmap_writable() as mm_geo:
...     for part_id in geofile.get_part_ids():
...         part = geofile.get_part_by_id(part_id)
...         arr = part.read_nodes(mm_geo)
...         arr[:, 0] += 1.0  # increment X coordiante

Creating geometry file from scratch

Please see tests/test_write_geometry.py, essentially you will need to do:

>>> from ensightreader import EnsightGeometryFile, GeometryPart, GeometryPart

>>> with open(output_geofile_path, "wb") as fp:
...     EnsightGeometryFile.write_header(fp)
...     GeometryPart.write_part_header(fp, part_id=1, part_name="TestElementTypes", node_coordinates=node_coordinates)
...     UnstructuredElementBlock.write_element_block(fp, element_type=et, connectivity=connectivity)

Code examples

The library comes with scripts that use it to convert EnSight Gold into other data formats. These are mostly useful for reference as they output data in text format, making them unsuitable for production use with large models.

ensight2obj script

This script converts surface elements of EnSight Gold parts into OBJ format (text). EnSight parts are represented as OBJ groups.

Demonstrates reading steady-state geometry, node coordinates, connectivity using traditional I/O.

For commandline usage, run the script with --help.

ensight2obj.ensight2obj(ensight_case_path: str, output_obj_path: str, part_name_regex: Optional[str] = None) int[source]

Main function of ensight2obj.py

ensight2vtk script

This script converts parts from EnSight Gold case into files in VTK legacy ASCII format.

Demonstrates reading steady-state geometry, node coordinates, connectivity, per-node and per-element variable data using memory-mapped I/O.

For commandline usage, run the script with --help.

ensight2vtk.ensight2vtk(ensight_case_path: str, output_vtk_path_given: str, part_name_regex: Optional[str] = None) int[source]

Main function of ensight2vtk.py

ensight2vtk.write_vtk_part(case: EnsightCaseFile, part_id: int, vtk_output_path: str) None[source]

Write part from EnSight Gold case with given ID as VTK legacy ASCII format file

ensight_transform script

This script does in-place transformation of node coordinates in given EnSight Gold case. Your original geofile will be modified!

Examples:

# increment X coordinate
ensight_transform --translate 1 0 0 sphere.case

# scale by 1000 (eg. m -> mm conversion)
ensight_transform --scale 1e3 1e3 1e3 sphere.case

# rotation matrix
ensight_transform --matrix \
    0 -1  0  0 \
    1  0  0  0 \
    0  0  1  0 \
    0  0  0  1 \
    sphere.case

# transform only "internalMesh" part
ensight_transform --translate 1 0 0 --only-parts internalMesh motorbike.case

For commandline usage, run the script with --help.

ensight_transform.ensight_transform(ensight_case_path: str, translate: Optional[ndarray[Any, dtype[float32]]] = None, scale: Optional[ndarray[Any, dtype[float32]]] = None, matrix: Optional[ndarray[Any, dtype[float32]]] = None, part_name_regex: Optional[str] = None) int[source]

Main function of ensight_transform.py