Utilities API
This document provides a reference to the osiris_utils utilities API, which includes tools for reading various OSIRIS data formats (Grid, Raw Particles, Tracks, HIST) and helper functions for common physics calculations and data operations.
Data Readers Structures
The package provides several classes for handling different types of OSIRIS data.
OsirisData
- class osiris_utils.data.data.OsirisData(filename)Source
Base class for handling OSIRIS simulation data files (HDF5, HIST and TIMINGS formats).
This class provides common functionality for reading and managing basic attributes from OSIRIS output files. It serves as the parent class for specialized data handlers.
- Parameters:
filename (str) – Path to the data file. Supported formats: - HDF5 files (.h5 extension) - HIST files (ending with _ene) - TIMINGS files (starting with timings)
- dtSource
Time step of the simulation [simulation units]
- Type:
float
- dimSource
Number of dimensions in the simulation (1, 2, or 3)
- Type:
int
- timeSource
Current simulation time and units as [value, unit_string]
- Type:
list[float, str]
- iterSource
Current iteration number
- Type:
int
- nameSource
Name identifier of the data field
- Type:
str
- typeSource
Type of data (e.g., ‘grid’, ‘particles’)
- Type:
str
- verboseSource
Verbosity flag controlling diagnostic messages (default: False)
- Type:
bool
Base class for all OSIRIS data types. All other data classes inherit from this class.
Key Attributes:
dt- Time step of the simulationdim- Dimensionality of the datatime- Physical time of the dataiter- Iteration numbername- Name of the datasettype- Type of the dataverbose- Verbosity flag for logging
- verbose(verbose: bool = True)Source
Set the verbosity of the class
- Parameters:
verbose (bool, optional) – If True, the class will print messages, by default True when calling (False when not calling)
Grid-Based Data
- class osiris_utils.data.data.OsirisGridFile(filename, data_slice: slice | None = None, load_data: bool = True)Source
Handles structured grid data from OSIRIS HDF5 simulations, including electromagnetic fields.
- Parameters:
filename (str) – Path to OSIRIS HDF5 grid file (.h5 extension)
- axisSource
Axis metadata with keys: - ‘name’: Axis identifier (e.g., ‘x1’) - ‘units’: Physical units (LaTeX formatted) - ‘long_name’: Descriptive name (LaTeX formatted) - ‘type’: Axis type (e.g., ‘SPATIAL’) - ‘plot_label’: Combined label for plotting
- Type:
list[dict]
Specialized class for handling grid-based field data such as electromagnetic fields.
Key Attributes:
grid- Grid informationnx- Number of grid pointsdx- Grid spacingx- Grid coordinatesaxis- Coordinate labelsdata- The actual field dataunits- Physical units of the datalabel- Data labels for visualization
Particle Data
- class osiris_utils.data.data.OsirisRawFile(filename)Source
Class to read the raw data from an OSIRIS HDF5 file.
- Parameters:
filename (str) – Path to OSIRIS HDF5 track file (.h5 extension)
Attributes
-----------
axis (dict[str, dict[str, str]]) –
- Dictionary where each key is a dataset name, and each value is another dictionary containing:
’name’ (str): Short name of the quantity (e.g., ‘x1’, ‘ene’)
’units’ (str): Units (LaTeX formatted, e.g., ‘c/omega_p’, ‘m_e c^2’)
’long_name’ (str): Descriptive name (LaTeX formatted, e.g., ‘x_1’, ‘En2’)
data (dict[str, np.ndarray]) – Dataset values indexed by dataset name (quants).
dim (int) – Number of spatial dimensions.
dt (float) – Time step between iterations.
grid (np.ndarray) – Grid boundaries as ((x1_min, x1_max), (x2_min, x2_max), …)
iter (int) – Iteration number corresponding to the data.
name (str) – Name of the species.
time (list[float, str]) – Simulation time and its units (e.g., [12.5, ‘1/omega_p’]).
type (str) – Type of data (e.g., ‘particles’ for raw files).
labels (list[str]) – Field labels/names (LaTeX formatted, e.g., ‘x_1’)
quants (list[str]) – field names of the data
units (list[str]) – Units of each field of the data (LaTeX formatted, e.g., ‘c/omega_p’)
Example
>>> import osiris_utils as ou >>> raw = ou.raw = ou.OsirisRawFile("path/to/raw/file.h5") >>> print(raw.data.keys()) >>> # Access x1 position of first 10 particles >>> print(raw.data["x1"][0:10]) >>> # Write beautiful labels and units >>> print("${} = $".format(raw.labels["x1"]) + "$[{}]$".format(track.units["x1"]))
Description: This class is responsible for reading and structuring raw particle data from an OSIRIS HDF5 file. It provides direct access to raw simulation data, in a structured dictionary format.
Key Attributes:
axis- Dictionary where each key is a dataset name, and each value is another dictionary containingname (str): The name of the quantity (e.g., r’x1’, r’ene’). units (str): The units associated with that dataset in LaTeX (e.g., r’c/\omega_p’, r’m_e c^2’). long_name (str): The name of the quantity in LaTeX (e.g., r’x_1’, r’En2’). dictionary of dictionaries
data- Dictionary mapping dataset names to their corresponding NumPy arrays.dim- Number of spatial dimensions in the simulation.dt- Simulation time step.grid- Array specifying the min/max coordinates of the simulation box along each axis.iter- Iteration number of the simulation snapshot.name- Name of the particle species in the dataset.time- List containing the simulation time and its associated units.type- Type of the dataset (e.g., “particles” for raw particle data).
Example Usage:
import osiris_utils as ou raw = ou.OsirisRawFile("path/to/raw/file.h5") print(raw.data.keys()) print(raw.data["x1"][0:10]) # Access x1 position of first 10 particles
Methods:
raw_to_file_tags(filename, type="all", n_tags=10, mask=None)Converts raw particle data into a file_tags format, used for selecting particles in OSIRIS tracking diagnostics.
Parameters:
filename(str): Path to the output tag file.type(str, optional): Selection mode (“all” for all tags, “random” for a random subset).n_tags(int, optional): Number of tags to select when type=”random”. Default is 10.mask(np.ndarray, optional): Boolean mask to filter valid particle tags before selection.
Output:
Generates a file storing selected particle tags to initialize OSIRIS track diagnostic.
Example Usage:
raw = ou.OsirisRawFile("path/to/raw/file/.../.h5") # Selecting 5 random tags from particles with energy>5 mask = raw.data["ene"] > 5. raw.raw_to_file_tags("output.tag", type="random", n_tags=5, mask=mask)
- raw_to_file_tags(filename, type: Literal['all', 'random'] = 'all', n_tags=10, mask=None)Source
Function to write a file_tags file from raw data. this file is used to choose particles for the OSIRIS track diagnostic.
- Parameters:
filename (str) – Path to the output file where tags will be stored.
type ({'all', 'random'}, optional) – Selection mode for tags: - ‘all’: Includes all available tags. - ‘random’: Randomly selects n_tags tags.
n_tags (int, optional) – Number of tags to randomly select when type is ‘random’. Default is 10.
mask (np.ndarray, optional) – Boolean mask array applied to filter valid tags before selection.
- Returns:
A file_tags file with path “filename” to be used for the OSIRIS track diagnostic.
Notes
The first element of the tag of a particle that is already being tracked is negative, so we apply the absolute function when generating the file
HIST Data
TIMINGS Data
TRACK Data
- class osiris_utils.data.data.OsirisTrackFile(filename)Source
Handles structured track data from OSIRIS HDF5 simulations.
- Parameters:
filename (str) – Path to OSIRIS HDF5 track file (.h5 extension)
- dataSource
dtype = [(field_name, float) for field_name in field_names] A structured numpy array with the track data Accessed as data[particles, time_iters][quant]
- Type:
numpy.ndarray of shape (num_particles, num_time_iter),
- num_particlesSource
Number of particlest tracked, they are accessed from 0 to num_particles-1
- Type:
int
- num_time_itersSource
Number of time iteratis, they are accessed from 0 to num_time_iters-1
- Type:
int
Example
>>> import osiris_utils as ou >>> track = ou.OsirisTrackFile(path/to/track_file.h5) >>> print(track.data[0:10, :]["x1"]) # Access x1 position of first 10 particles over all time steps
Specialized class for handling particle track data from OSIRIS HDF5 simulations. This class processes and organizes particle tracking data into a structured numpy array, providing an intuitive interface for accessing particle properties over time.
Key Attributes:
grid- Grid boundaries of the simulationdata- Structured numpy array with tracked particle data, accessed asdata[particle, time_iter][quant]units- List of units corresponding to each tracked fieldlabels- LaTeX-formatted labels for each tracked quantityquants- Field names of the tracked datanum_particles- Number of tracked particles, indexed from0tonum_particles - 1num_time_iters- Number of time iterations, indexed from0tonum_time_iters - 1
Methods:
__array__()- Returns thedataattribute as a standard numpy array.__str__()- Provides a string summary of the object, including simulation grid, iteration, and dimensions._load_basic_attributes(f: h5py.File)- Loads essential attributes such as simulation time step (dt), dimensionality (dim), and file metadata.
Example Usage:
import osiris_utils as ou track = ou.OsirisTrackFile("path/to/track_file.h5") print(track.data[0:10, :]["x1"]) # Access x1 position of first 10 particles over all time steps
Convert track file to the older more readable format
- convert_tracks(filename_in)
Description: Converts a new OSIRIS track file aka IDL-formatted aka tracks-2 to an older format that is more human-readable. In the old format, each particle is stored in a separate folder, with datasets for each quantity. The function reads data from the input file, processes it, and writes it to a new file with the suffix “-v2”.
Parameters: -
filename_in(str): The path to the input IDL-formatted track file.Output File: - A new HDF5 file is created with the same name as the input file, but with the suffix -v2 added to the filename. This file will be in the old, more readable format.
Example:
>>> import osiris_utils as ou >>> ou.convert_tracks('path/to/input_trackfile.h5') >>> # The output will be saved as 'path/to/input_trackfile-v2.h5'
Notes: - The old format stores each particle’s data in separate groups within the file. - This function assumes the input file follows the IDL (new) format as expected by OSIRIS. - Code from https://github.com/GoLP-IST/RaDi-x/blob/main/tools/convert_idl_tracks_py3.py
Physics & Analysis
Numerical Methods
Methods for common physics calculations:
courant2D- Calculate the Courant condition for 2D simulationstransverse_average- Compute averages along transverse directions
- osiris_utils.utils.courant2D(dx: float, dy: float) floatSource
Compute the Courant number for a 2D simulation.
- Parameters:
dx (float) – The spacing in the x direction.
dy (float) – The spacing in the y direction.
- Returns:
float – The limit for dt.
- osiris_utils.utils.transverse_average(data: ndarray) ndarraySource
Computes the transverse average of a 2D array.
- Parameters:
data (numpy.ndarray) – Dim: 2D. The input data.
- Returns:
numpy.ndarray – Dim: 1D. The transverse average.
Data Utilities
File Operations
Utilities for data handling and file operations:
time_estimation- Estimate runtime for operationsfilesize_estimation- Estimate file sizesintegrate- Numerical integration routinessave_data- Save data to diskread_data- Read data from disk
- osiris_utils.utils.time_estimation(n_cells: int, ppc: int, t_steps: int, n_cpu: int, push_time: float = 1e-07, hours: bool = False) floatSource
Estimate the simulation time.
- Parameters:
n_cells (int) – The number of cells.
ppc (int) – The number of particles per cell.
push_time (float) – The time per push.
t_steps (int) – The number of time steps.
n_cpu (int) – The number of CPU’s.
hours (bool, optional) – If True, the output will be in hours. If False, the output will be in seconds. The default is False.
- Returns:
float – The estimated time in seconds or hours.
- osiris_utils.utils.integrate(array: ndarray, dx: float, axis: int = None, start_side: str = 'left') ndarraySource
Integrate a N-D array using the cumulative Simpson’s rule, from right to left along a given axis
- Parameters:
array (numpy.ndarray) – The input array.
dx (float) – The spacing between points.
axis (int, optional) – The axis along which to integrate. If None, the default is 0. This will assume that the input array is 1D.
- osiris_utils.utils.save_data(data: ndarray, savename: str, option: str = 'numpy') NoneSource
Save the data to a .txt (with Numpy) or .csv (with Pandas) file.
- Parameters:
data (list) – The data to be saved.
savename (str) – The path to the file.
option (str, optional) – The option for saving the data. The default is ‘numpy’. Can be ‘numpy’ or ‘pandas’.
- osiris_utils.utils.read_data(filename: str, option: str = 'numpy') ndarraySource
Read the data from a .txt file.
- Parameters:
filename (str) – The path to the file.
- Returns:
numpy.ndarray – Dim: 2D. The data.