mhd_io#

Lazy HDF readers for PSI MAS and POT3D magnetohydrodynamic model output.

This module provides a unified, unit-aware interface for reading three-dimensional MHD model output from Predictive Science Inc.’s MAS and POT3D solvers. Both HDF4 (.hdf) and HDF5 (.h5) files are supported through a common API.

Class hierarchy#

The implementation uses a mixin-based hierarchy that separates file format concerns from physical model concerns:

_HdfInterface (ABC)                 ← public interface contract
└── _HdfData (_HdfInterface, ABC)   ← shared MAS/POT3D logic
    ├── H5MasData  (_H5DataMixin, _HdfData)   ← HDF5 MAS output
    ├── H4MasData  (_H4DataMixin, _HdfData)   ← HDF4 MAS output
    ├── H5Pot3dData(_H5DataMixin, _HdfData)   ← HDF5 POT3D output
    └── H4Pot3dData(_H4DataMixin, _HdfData)   ← HDF4 POT3D output

_HdfInterface (ABC)
└── _HdfScale (_HdfInterface, ABC)  ← coordinate axis arrays
    ├── H5Scale                     ← HDF5 dimension scale
    └── H4Scale                     ← HDF4 SDS dimension

Format mixins inject all file-I/O and raw-array logic:

_H5DataMixin ← h5py-backed I/O (open/close/data property) _H4DataMixin ← pyhdf-backed I/O (open/close/data property)

Lazy-loading pattern#

All concrete data classes open the HDF file and read metadata at instantiation, but do not load the array into memory until _HdfData.read() is called. This makes it cheap to inspect many files:

>>> from pathlib import Path
>>> from psi_io.mhd_io import PsiData
>>> reader = PsiData('br001001.h5', model='mas') # file opened, metadata parsed
>>> reader.quantity                              # 'br' — from filename / file attrs
>>> data, r, t, p = reader.read()               # array loaded here

Axis ordering convention#

PSI HDF files are written by Fortran code in column-major order, so when h5py reads them into numpy (row-major), the physical (r, θ, φ) coordinate ordering is reversed to (φ, θ, r) in the numpy array. The reader hides this inversion:

  • _HdfData.shape reports the numpy storage shape (Nφ, Nθ, Nr).

  • _HdfInterface.__getitem__() accepts indices in the physical (r, θ, φ) user order and internally reverses them to match numpy storage.

Public factory#

PsiData() is the recommended entry point. It inspects the file extension and the model keyword to select the correct concrete class:

>>> reader = PsiData('br001001.h5', model='mas')
>>> reader = PsiData('br001001.hdf', model='pot3d')

See Also#

psi_io._props : Physical property metadata (units, mesh stagger) for each quantity. psi_io._units : Physical normalization constants and custom astropy units. psi_io._mesh : Staggered-grid utilities used during mesh conversion on read.

Functions

extract_quantity_from_filepath(ifile[, default])

Extract the MAS/POT3D quantity name from a filename stem.

extract_sequence_from_filepath(ifile[, default])

Extract the sequence number from a filename stem.

get_mas_quantity_properties(variable)

Return the Props descriptor for a MAS quantity.

get_pot3d_quantity_properties(variable)

Return the Props descriptor for a POT3D quantity.

get_psi_scale_properties(variable)

Return the Props descriptor for a PSI coordinate scale.

PsiData(ifile, /[, model])

Open a PSI MAS or POT3D HDF file and return the appropriate data reader.

Attributes

ModelType

Literal type alias for the three recognized PSI model types.