API Reference
Case file
- ensightreader.read_case(path: str) EnsightCaseFile [source]
Read EnSight Gold case
This will only parse the case file, not any data files. Use the returned object to load whichever data you need.
- Parameters
path – Path to the
*.case
file- Returns
EnsightCaseFile
object
- class ensightreader.EnsightCaseFile(casefile_path: str, geometry_model: ~ensightreader.EnsightGeometryFileSet, variables: ~typing.Dict[str, ~ensightreader.EnsightVariableFileSet], constant_variables: ~typing.Dict[str, ~ensightreader.EnsightConstantVariable], timesets: ~typing.Dict[int, ~ensightreader.Timeset], _geometry_file_cache: ~typing.Dict[int, ~ensightreader.EnsightGeometryFile] = <factory>, _variable_file_cache: ~typing.Dict[~typing.Tuple[int, str], ~ensightreader.EnsightVariableFile] = <factory>)[source]
EnSight Gold case file
To load a case, use:
>>> import ensightreader >>> case = ensightreader.read_case("example.case") >>> geofile = case.get_geometry_model() >>> velocity_variable = case.get_variable("U")
EnSight Gold casefile is text file with description of the case and links to data files where the actual geometry and variables are stored.
Since the case may be transient and therefore have multiple data files for the same variable, etc., there is an indirection: the actual data is accessed using
EnsightGeometryFile
andEnsightVariableFile
objects, which are constructed usingEnsightGeometryFileSet
,EnsightVariableFileSet
helper objects insideEnsightCaseFile
.Note
for geometry, only
model
is supported (nomeasured
,match
orboundary
)only unstructured grid (
coordinates
) geometry is supported; structured parts are not supportedfilesets and single-file cases are not supported; only
C Binary
files are supportedchange_coords_only
is not supportedchanging_geometry_per_part
is only read and stored, not handled in reading geometry datasome variable types are unsupported, see
VariableType
ghost cells are not supported
cases with more than one time set are supported, but you may not be able to use
EnsightCaseFile.get_variable()
since the library expects geometry and variable time sets to match; if you know which geofile belongs to which variable file you can read the variable files usingEnsightVariableFileSet.get_file()
.
- casefile_path
path to the
*.case
file- Type
str
- geometry_model
accessor to the
model
geometry data
- variables
dictionary mapping variable names keys to variable data accessors
- Type
Dict[str, ensightreader.EnsightVariableFileSet]
- constant_variables
dictionary mapping constant variable names keys to
EnsightConstantVariable
- Type
Dict[str, ensightreader.EnsightConstantVariable]
- timesets
dictionary mapping time set IDs to time set description objects
- Type
Dict[int, ensightreader.Timeset]
- classmethod from_file(casefile_path: str) EnsightCaseFile [source]
Read EnSight Gold case
Note
This only reads the casefile, not any data files. Any geometry or variable data that you want to read must be read explicitly.
- Parameters
casefile_path – path to the
*.case
file
- get_constant_variable_value(name: str, timestep: int = 0) float [source]
Get value of constant per-case variable for given timestep
- Parameters
name – name of the constant variable
timestep – number of timestep, starting from zero (for non-transient variable, use 0)
- get_geometry_model(timestep: int = 0) EnsightGeometryFile [source]
Get geometry for given timestep
Note
The returned
EnsightGeometryFile
is cached, so it will not parse the file again when you request the geometry again.- Parameters
timestep – number of timestep, starting from zero (for non-transient geometry, use 0)
- get_time_values(timeset_id: Optional[int] = None) Optional[List[float]] [source]
Return time values for given time set
This will return time values as specified in the case file (this could be physical time in seconds, or something else). Timesteps
0
,1
, … occur at timescase.get_time_values()[0]
,case.get_time_values()[1]
, …- Parameters
timeset_id – time set ID for which to return time values - if there are less than two time sets, you can omit this parameter
- Returns
List of time values, or None if there are no time sets defined
- get_variable(name: str, timestep: int = 0) EnsightVariableFile [source]
Get variable for given timestep
Note
The returned
EnsightVariableFile
is cached, so it will not parse the file again when you request the variable again.- Parameters
name – name of the variable
timestep – number of timestep, starting from zero (for non-transient variable, use 0)
- get_variables() List[str] [source]
Return list of variable names
Note that this does not include per-case constants.
- is_transient() bool [source]
Return True if the case is transient (has at least one time set defined)
- to_file(casefile_path: str) None [source]
Write EnSight Gold case to file
See
EnsightCaseFile.to_string()
for more.- Parameters
casefile_path – path to the
*.case
file
- to_string() str [source]
Return .case file contents as a string
This method is useful if you wish to modify the case, eg. add/remove variables, apply offset to time values, etc.
Note
This method works by serializing the internal representation, meaning that any lines in the original .case file that were skipped due to missing support from the library (as well as comments) will not appear in the output at all.
- class ensightreader.EnsightGeometryFileSet(casefile_dir_path: str, timeset: Optional[Timeset], filename: str, changing_geometry_per_part: bool = False)[source]
Helper object for loading geometry files
This is used by
EnsightCaseFile.get_geometry_model()
to give youEnsightGeometryFile
.- casefile_dir_path
path to casefile directory (root for relative paths in casefile)
- Type
str
- timeset
time set in which the geometry is defined, or None if it’s not transient
- Type
Optional[ensightreader.Timeset]
- filename
path to the data file(s), including
*
wildcards if transient- Type
str
- changing_geometry_per_part
whether parts contain information about type of transient changes
- Type
bool
- get_file(timestep: int = 0) EnsightGeometryFile [source]
Return geometry for given timestep (use 0 if not transient)
- class ensightreader.EnsightVariableFileSet(casefile_dir_path: str, timeset: Optional[Timeset], variable_location: VariableLocation, variable_type: VariableType, variable_name: str, filename: str)[source]
Helper object for loading variable files
This is used by
EnsightCaseFile.get_variable()
to give youEnsightVariableFile
.- casefile_dir_path
path to casefile directory (root for relative paths in casefile)
- Type
str
- timeset
time set in which the variable is defined, or None if it’s not transient
- Type
Optional[ensightreader.Timeset]
- variable_location
where the variable is defined (elements or nodes)
- variable_type
type of the variable (scalar, …)
- variable_name
name of the variable (‘description’ field in casefile)
- Type
str
- filename
path to the data file(s), including
*
wildcards if transient- Type
str
- get_file(geofile: EnsightGeometryFile, timestep: int = 0) EnsightVariableFile [source]
Return variable for given timestep (use 0 if not transient)
Note
Due to how EnSight Gold format works, you need to have matching geofile already parsed before you attempt to read the variable file. The variable file itself (alas) does not give lengths of the arrays it contains.
For transient geometry combined with transient variable, pay extra care because they have to match - different timesteps can have different number of elements, nodes, etc.
- class ensightreader.Timeset(timeset_id: int, description: Optional[str], number_of_steps: int, filename_numbers: List[int], time_values: List[float])[source]
Description of time set in EnSight Gold case
This means a non-decreasing sequence of times for which geometry and/or variable values are known.
- timeset_id
ID of the time set
- Type
int
- description
label of the time set, or None
- Type
Optional[str]
- number_of_steps
number of timesteps
- Type
int
- filename_numbers
list of numbers for filenames (to be filled in place of
*
wildcards)- Type
List[int]
- time_values
list of time values (ie. seconds, or something else)
- Type
List[float]
Geometry file
- class ensightreader.EnsightGeometryFile(file_path: str, description_line1: str, description_line2: str, node_id_handling: IdHandling, element_id_handling: IdHandling, extents: Optional[ndarray[Any, dtype[float32]]], parts: Dict[int, GeometryPart], changing_geometry_per_part: bool)[source]
EnSight Gold binary geometry file
To use it:
>>> import ensightreader >>> case = ensightreader.read_case("example.case") >>> geofile = case.get_geometry_model() >>> part_names = geofile.get_part_names() >>> part = geofile.get_part_by_name(part_names[0])
This objects contains metadata parsed from EnSight Gold geometry file; it does not hold any node coordinates, connectivity, etc., but after it’s been created it can read any such data on demand.
- file_path
path to the actual geometry file (no wildcards)
- Type
str
- description_line1
first line in header
- Type
str
- description_line2
second line in header
- Type
str
- node_id_handling
node ID presence
- element_id_handling
element ID presence
- extents
optional extents given in header (xmin, xmax, ymin, ymax, zmin, zmax)
- Type
Optional[numpy.ndarray[Any, numpy.dtype[numpy.float32]]]
- parts
dictonary mapping part IDs to
GeometryPart
objects- Type
Dict[int, ensightreader.GeometryPart]
- changing_geometry_per_part
whether parts contain information about type of transient changes
- Type
bool
- classmethod from_file_path(file_path: str, changing_geometry_per_part: bool) EnsightGeometryFile [source]
Parse EnSight Gold geometry file
- get_part_by_id(part_id: int) Optional[GeometryPart] [source]
Return part with given part ID, or None
- get_part_by_name(name: str) Optional[GeometryPart] [source]
Return part with given name, or None
- mmap() Generator[mmap, None, None] [source]
Return read-only memorymap of the file (convenience method)
Use this in
with
block and pass the resulting object toGeometryPart
methods.Note
This is preferred over the similar
open()
method if you don’t need a copy of the data since the returned arrays will be backed by the memorymap. Be careful to keep the memorymap around as long as you need the arrays.Usage:
with geometry_file.mmap() as fp_geo: nodes = part.read_nodes(fp_geo)
Equivalent code:
with open(geometry_file.file_path, "rb") as fp_geo, _mmap.mmap(fp_geo.fileno(), 0, access=mmap.ACCESS_READ) as mm_geo: nodes = part.read_nodes(mm_geo)
- mmap_writable() Generator[mmap, None, None] [source]
Return writable memorymap of the file (convenience method)
Use this in
with
block and pass the resulting object toGeometryPart
methods.Note
This special version of the
mmap()
method can be used if you wish to modify the underlying file. Use carefully.Usage:
with geometry_file.mmap_writable() as fp_geo: nodes = part.read_nodes(fp_geo) nodes[:, 0] = 0.0 # set first coordinate to zero for part nodes
Equivalent code:
with open(geometry_file.file_path, "r+b") as fp_geo, _mmap.mmap(fp_geo.fileno(), 0, access=mmap.ACCESS_WRITE) as mm_geo: nodes = part.read_nodes(mm_geo) nodes[:, 0] = 0.0 # set X coordinate to zero for part nodes
- open() BinaryIO [source]
Return the opened file in read-only binary mode (convenience method)
Use this in
with
block and pass the resulting object toGeometryPart
methods.Note
This is a simpler alternative to the similar
mmap()
method - it doesn’t require you to reason about lifetime of the opened file since you will be getting arrays which own their data (as opposed to views into the memory-mapped file).Usage:
with geometry_file.open() as fp_geo: nodes = part.read_nodes(fp_geo)
Equivalent code:
with open(geometry_file.file_path, "rb") as fp_geo: nodes = part.read_nodes(fp_geo)
- static write_header(fp: Union[BinaryIO, mmap], description_line1: str = 'Generated by ensightreader', description_line2: str = '', node_id_handling: IdHandling = IdHandling.OFF, element_id_handling: IdHandling = IdHandling.OFF, extents: Optional[ndarray[Any, dtype[float32]]] = None) None [source]
Writes geometry file header to given opened file
Make sure that
fp
seek position is at the beginning of the file.
- class ensightreader.GeometryPart(offset: int, part_id: int, part_name: str, number_of_nodes: int, element_blocks: List[UnstructuredElementBlock], node_id_handling: IdHandling, element_id_handling: IdHandling, changing_geometry: Optional[ChangingGeometry] = None)[source]
A part in EnSight Gold geometry file
To use it:
>>> import ensightreader >>> case = ensightreader.read_case("example.case") >>> geofile = case.get_geometry_model() >>> part_names = geofile.get_part_names() >>> part = geofile.get_part_by_name(part_names[0]) >>> print(part.is_volume()) >>> print(part.number_of_nodes)
Each geometry part has its own set of nodes not shared with other parts. Elements are defined in blocks, where each block may have different type (tetrahedra, wedge, etc.) but all elements in the block have the same type.
- offset
offset to ‘part’ line in the geometry file
- Type
int
- part_id
part number
- Type
int
- part_name
part name (‘description’ line)
- Type
str
- number_of_nodes
number of nodes for this part
- Type
int
- element_blocks
list of element blocks, in order of definition in the file
- Type
- node_id_handling
node ID presence
- element_id_handling
element ID presence
- changing_geometry
type of transient changes (coordinate/connectivity/none)
- Type
Optional[ensightreader.ChangingGeometry]
- classmethod from_file(fp: Union[BinaryIO, mmap], node_id_handling: IdHandling, element_id_handling: IdHandling, changing_geometry_per_part: bool) GeometryPart [source]
Used internally by
EnsightGeometryFile.from_file_path()
- get_number_of_elements_of_type(element_type: ElementType) int [source]
Return number of elements (of given type)
- property number_of_elements: int
Return number of elements (of all types)
- read_node_ids(fp: Union[BinaryIO, mmap]) Optional[ndarray[Any, dtype[int32]]] [source]
Read node IDs for this part, if present
Note
This method simply returns the node ID array as present in the file; it does not differentiate between
node id given
andnode id ignore
, etc.- Parameters
fp – opened geometry file object in
"rb"
mode- Returns
1D array of int32 with node IDs, or None if node IDs are not present in the file
- read_nodes(fp: Union[BinaryIO, mmap]) ndarray[Any, dtype[float32]] [source]
Read node coordinates for this part
- Parameters
fp – opened geometry file object in
"rb"
mode- Returns
2D
(n, 3)
array of float32 with node coordinates
- static write_part_header(fp: Union[BinaryIO, mmap], part_id: int, part_name: str, node_coordinates: ndarray[Any, dtype[float32]], node_ids: Optional[ndarray[Any, dtype[int32]]] = None) None [source]
Write part header to given opened file
node_coordiantes
are supposed to be (N, 3)-shaped float32 array.
- class ensightreader.UnstructuredElementBlock(offset: int, number_of_elements: int, element_type: ElementType, element_id_handling: IdHandling, part_id: int)[source]
A block of elements of the same type in a part in EnSight Gold binary geometry file
To use it:
>>> from ensightreader import read_case, ElementType >>> case = read_case("example.case") >>> geofile = case.get_geometry_model() >>> part_names = geofile.get_part_names() >>> part = geofile.get_part_by_name(part_names[0]) >>> with open(geofile.file_path, "rb") as fp_geo: ... for block in part.element_blocks: ... if block.element_type == ElementType.NFACED: ... polyhedra_face_counts, face_node_counts, face_connectivity = block.read_connectivity_nfaced(fp_geo) ... elif block.element_type == ElementType.NSIDED: ... polygon_node_counts, polygon_connectivity = block.read_connectivity_nsided(fp_geo) ... else: ... connectivity = block.read_connectivity(fp_geo)
- offset
offset to ‘element type’ line in file (eg. ‘tria3’)
- Type
int
- number_of_elements
number of elements in this block
- Type
int
- element_type
type of elements in this block
- element_id_handling
element ID presence
- part_id
part number
- Type
int
- classmethod from_file(fp: Union[BinaryIO, mmap], element_id_handling: IdHandling, part_id: int) UnstructuredElementBlock [source]
Used internally by
GeometryPart.from_file()
- read_connectivity(fp: Union[BinaryIO, mmap]) ndarray[Any, dtype[int32]] [source]
Read connectivity (for elements other than NSIDED/NFACED)
Use this for elements which have constant number of nodes per element (ie. any element type except polygons and polyhedra).
- Parameters
fp – opened geometry file object in
"rb"
mode- Returns
2D
(n, k)
array of int32 with node indices (numbered from 1), wheren
is the number of elements andk
is the number of nodes defining each element
- read_connectivity_nfaced(fp: Union[BinaryIO, mmap]) Tuple[ndarray[Any, dtype[int32]], ndarray[Any, dtype[int32]], ndarray[Any, dtype[int32]]] [source]
Read connectivity (for NFACED elements)
- Parameters
fp – opened geometry file object in
"rb"
mode- Returns
tuple
(polyhedra_face_counts, face_node_counts, face_connectivity)
wherepolyhedra_face_counts
is 1D array of type int32 giving number of faces for each polygon,face_node_counts
is 1D array of type int32 giving number of nodes for each face for each polygon (ordered as 1st polygon 1st face, 1st polygon 2nd face, …, 2nd polygon 1st face, …),face_connectivity
is 1D array of type int32 giving node indices (numbered from 1)
- read_connectivity_nsided(fp: Union[BinaryIO, mmap]) Tuple[ndarray[Any, dtype[int32]], ndarray[Any, dtype[int32]]] [source]
Read connectivity (for NSIDED elements)
- Parameters
fp – opened geometry file object in
"rb"
mode- Returns
tuple
(polygon_node_counts, polygon_connectivity)
wherepolygon_node_counts
is 1D array of type int32 giving number of nodes for each polygon andpolygon_connectivity
is 1D array of type int32 giving node indices (numbered from 1) for every polygon
- read_element_ids(fp: Union[BinaryIO, mmap]) Optional[ndarray[Any, dtype[int32]]] [source]
Read element IDs for this element block, if present
Note
This method simply returns the element ID array as present in the file; it does not differentiate between
element id given
andelement id ignore
, etc.- Parameters
fp – opened geometry file object in
"rb"
mode- Returns
1D array of int32 with element IDs, or None if element IDs are not present in the file
- static write_element_block(fp: Union[BinaryIO, mmap], element_type: ElementType, connectivity: ndarray[Any, dtype[int32]], element_ids: Optional[ndarray[Any, dtype[int32]]] = None) None [source]
Write element block (not NSIDED/NFACED) to given opened file
- static write_element_block_nfaced(fp: Union[BinaryIO, mmap], polyhedra_face_counts: ndarray[Any, dtype[int32]], face_node_counts: ndarray[Any, dtype[int32]], face_connectivity: ndarray[Any, dtype[int32]], element_ids: Optional[ndarray[Any, dtype[int32]]] = None) None [source]
Write NFACED element block to given opened file
Variable file
- class ensightreader.EnsightVariableFile(file_path: str, description_line: str, variable_name: str, variable_location: VariableLocation, variable_type: VariableType, part_offsets: Dict[int, int], part_element_offsets: Dict[Tuple[int, ElementType], int], geometry_file: EnsightGeometryFile, part_per_node_undefined_values: Dict[int, float], part_per_element_undefined_values: Dict[Tuple[int, ElementType], float])[source]
EnSight Gold binary variable file
To use it:
>>> import ensightreader >>> case = ensightreader.read_case("example.case") >>> geofile = case.get_geometry_model() >>> velocity_variable = case.get_variable("U") >>> part_names = geofile.get_part_names() >>> part = geofile.get_part_by_name(part_names[0]) >>> with open(velocity_variable.file_path, "rb") as fp_var: ... part_velocity = velocity_variable.read_node_data(fp_var, part.part_id)
Note
there are some limitations for per-element variable files, see
EnsightVariableFile.read_element_data()
coordinates partial
is not supported
This objects contains metadata parsed from EnSight Gold variable file; it does not hold any variable data arrays, but after it’s been created it can read any such data on demand.
- file_path
path to the actual variable file (no wildcards)
- Type
str
- description_line
line in header
- Type
str
- variable_name
name of the variable in case file
- Type
str
- variable_location
where the variable is defined (elements or nodes)
- variable_type
type of the variable (scalar, …)
- part_offsets
dictionary mapping part IDs to offset to ‘part’ line in file
- Type
Dict[int, int]
- part_element_offsets
for per-element variables, this holds a dictionary mapping
(part ID, element type)
tuples to offset to ‘element type’ line in file- Type
Dict[Tuple[int, ensightreader.ElementType], int]
- part_per_node_undefined_values
for per-node variables, this holds a dictionary mapping part IDs to value that should be considered as undefined (
coordinates undef
)- Type
Dict[int, float]
- part_per_element_undefined_values
for per-element variables, this holds a dictionary mapping
(part ID, element type)
tuples to value that should be considered as undefined (element_type undef
)- Type
Dict[Tuple[int, ensightreader.ElementType], float]
- classmethod from_file_path(file_path: str, variable_name: str, variable_location: VariableLocation, variable_type: VariableType, geofile: EnsightGeometryFile) EnsightVariableFile [source]
Used internally by
EnsightVariableFileSet.get_file()
- is_defined_for_part_id(part_id: int) bool [source]
Return True if variable is defined for given part, else False
- mmap() Generator[mmap, None, None] [source]
Return read-only memorymap of the file (convenience method)
Use this in
with
block and pass the resulting object toread_node_data()
andread_element_data()
methods.Note
This is preferred over the similar
open()
method if you don’t need a copy of the data since the returned arrays will be backed by the memorymap. Be careful to keep the memorymap around as long as you need the arrays.Usage:
with variable_file.mmap() as mm_var: variable_data = variable.read_node_data(mm_var, part.part_id)
Equivalent code:
with open(variable_file.file_path, "rb") as fp_var, _mmap.mmap(fp_var.fileno(), 0, access=mmap.ACCESS_READ) as mm_var: variable_data = variable.read_node_data(mm_var, part.part_id)
- mmap_writable() Generator[mmap, None, None] [source]
Return writable memorymap of the file (convenience method)
Use this in
with
block and pass the resulting object toread_node_data()
andread_element_data()
methods.Note
This special version of the
mmap()
method can be used if you wish to modify the underlying file. Use carefully.Usage:
with geometry_file.mmap_writable() as fp_geo: nodes = part.read_nodes(fp_geo) nodes[:, 0] = 0.0 # set X coordinate to zero for part nodes
Equivalent code:
with open(geometry_file.file_path, "r+b") as fp_geo, _mmap.mmap(fp_geo.fileno(), 0, access=mmap.ACCESS_WRITE) as mm_geo: nodes = part.read_nodes(mm_geo)
Usage:
with variable_file.mmap_writable() as mm_var: variable_data = variable.read_node_data(mm_var, part.part_id) variable_data[:] = np.sqrt(variable_data) # apply square root function to the data
Equivalent code:
with open(variable_file.file_path, "r+b") as fp_var, _mmap.mmap(fp_var.fileno(), 0, access=mmap.ACCESS_WRITE) as mm_var: variable_data = variable.read_node_data(mm_var, part.part_id) variable_data[:] = np.sqrt(variable_data) # apply square root function to the data
- open() BinaryIO [source]
Return the opened file in read-only binary mode (convenience method)
Use this in
with
block and pass the resulting object toread_node_data()
andread_element_data()
methods.Note
This is a simpler alternative to the similar
mmap()
method - it doesn’t require you to reason about lifetime of the opened file since you will be getting arrays which own their data (as opposed to views into the memory-mapped file).Usage:
with variable_file.open() as fp_var: variable_data = variable.read_node_data(fp_var, part.part_id)
Equivalent code:
with open(variable_file.file_path, "rb") as fp_var: variable_data = variable.read_node_data(fp_var, part.part_id)
- read_element_data(fp: Union[BinaryIO, mmap], part_id: int, element_type: ElementType) Optional[ndarray[Any, dtype[float32]]] [source]
Read per-element variable data for given part and element type
Due to how EnSight Gold format works, variable file mirrors the geometry file in that per-element values are (alas) defined separately for each element block. This introduces edge cases (missing element blocks so that the variable could only be partially defined for the part; multiple element blocks of the same element type; etc.) which may be tedious to handle. This implementation relies on the assumption that there are no repeated blocks of the same type.
Note
Variable is always either per-node or per-element; be sure to call the right method for your data.
- Parameters
fp – opened variable file object in “rb” mode
part_id – part number for which to read data
element_type – element type for which to read data (typically, you want to iterate over element blocks in the part and retrieve data for their respective element types)
- Returns
If the variable is not defined for the part, None is returned. If the variable is defined and is a scalar, 1D array of float32 is returned. Otherwise the returned value is 2D
(n, k)
array of float32 wheren
is number of nodes andk
is number of values based on variable type (vector, tensor).
- read_node_data(fp: Union[BinaryIO, mmap], part_id: int) Optional[ndarray[Any, dtype[float32]]] [source]
Read per-node variable data for given part
Note
Variable is always either per-node or per-element; be sure to call the right method for your data.
- Parameters
fp – opened variable file object in “rb” mode
part_id – part number for which to read data
- Returns
If the variable is not defined for the part, None is returned. If the variable is defined and is a scalar, 1D array of float32 is returned. Otherwise the returned value is 2D
(n, k)
array of float32 wheren
is number of nodes andk
is number of values based on variable type (vector, tensor).
Enum types
- enum ensightreader.ElementType(value)[source]
Element type in EnSight Gold geometry file
Note
Ghost cell variants
g_*
are not supported.Valid values are as follows:
- POINT = <ElementType.POINT: 'point'>
- BAR2 = <ElementType.BAR2: 'bar2'>
- BAR3 = <ElementType.BAR3: 'bar3'>
- TRIA3 = <ElementType.TRIA3: 'tria3'>
- TRIA6 = <ElementType.TRIA6: 'tria6'>
- QUAD4 = <ElementType.QUAD4: 'quad4'>
- QUAD8 = <ElementType.QUAD8: 'quad8'>
- TETRA4 = <ElementType.TETRA4: 'tetra4'>
- TETRA10 = <ElementType.TETRA10: 'tetra10'>
- PYRAMID5 = <ElementType.PYRAMID5: 'pyramid5'>
- PYRAMID13 = <ElementType.PYRAMID13: 'pyramid13'>
- PENTA6 = <ElementType.PENTA6: 'penta6'>
- PENTA15 = <ElementType.PENTA15: 'penta15'>
- HEXA8 = <ElementType.HEXA8: 'hexa8'>
- HEXA20 = <ElementType.HEXA20: 'hexa20'>
- NSIDED = <ElementType.NSIDED: 'nsided'>
- NFACED = <ElementType.NFACED: 'nfaced'>
The
Enum
and its members also have the following methods:- property dimension: int
Return dimension of element
Returns 3 for volume elements, 2 for surface elements, 1 for line elements and 0 for point elements.
- property nodes_per_element: int
Return number nodes defining the element
This only makes sense for elements consisting of constant number of nodes. For NSIDED and NFACED element type, this raises and exception.
- enum ensightreader.VariableType(value)[source]
Type of variable in EnSight Gold case
Note
Complex variables and “per measured” variables are not supported.
Valid values are as follows:
- SCALAR = <VariableType.SCALAR: 'scalar'>
- VECTOR = <VariableType.VECTOR: 'vector'>
- TENSOR_SYMM = <VariableType.TENSOR_SYMM: 'tensor symm'>
- TENSOR_ASYM = <VariableType.TENSOR_ASYM: 'tensor asym'>
- enum ensightreader.VariableLocation(value)[source]
Location of variable in EnSight Gold case
Whether the variable is defined for cells or nodes.
Valid values are as follows:
- PER_ELEMENT = <VariableLocation.PER_ELEMENT: 'element'>
- PER_NODE = <VariableLocation.PER_NODE: 'node'>
- enum ensightreader.IdHandling(value)[source]
Handling of node/element IDs in EnSight Gold geometry file.
This is defined in geometry file header and describes whether IDs are present in the file or not.
Valid values are as follows:
- OFF = <IdHandling.OFF: 'off'>
- GIVEN = <IdHandling.GIVEN: 'given'>
- ASSIGN = <IdHandling.ASSIGN: 'assign'>
- IGNORE = <IdHandling.IGNORE: 'ignore'>
The
Enum
and its members also have the following methods:- property ids_present: bool
Return True if IDs are present in geometry file, otherwise False
- enum ensightreader.ChangingGeometry(value)[source]
Additional information about transient geometry
Valid values are as follows:
- NO_CHANGE = <ChangingGeometry.NO_CHANGE: 'no_change'>
- COORD_CHANGE = <ChangingGeometry.COORD_CHANGE: 'coord_change'>
- CONN_CHANGE = <ChangingGeometry.CONN_CHANGE: 'conn_change'>
Exceptions and warnings
- class ensightreader.EnsightReaderError(msg: str, fp: Optional[Union[TextIO, BinaryIO, mmap]] = None, lineno: Optional[int] = None)[source]
Error raised when parsing EnSight Gold binary files
- file_path
path to file where the error was encountered
- Type
str
- file_offset
approximate seek position of the error (this may be a bit past the place where the error is - it’s the seek position when this exception was raised)
- Type
int
- file_lineno
line number of the error (this only applies to errors in
*.case
file, as other files are binary)- Type
int
- class ensightreader.EnsightReaderWarning(msg: str, fp: Optional[Union[TextIO, BinaryIO, mmap]] = None, lineno: Optional[int] = None)[source]
Warning raised when parsing EnSight Gold binary files
- file_path
path to file where the error was encountered
- Type
str
- file_offset
approximate seek position of the error (this may be a bit past the place where the error is - it’s the seek position when this exception was raised)
- Type
int
- file_lineno
line number of the error (this only applies to errors in
*.case
file, as other files are binary)- Type
int