API Reference

Python module definition

taco_format.calc_rmsd(coords1, coords2)

Calculate root-mean-square deviation (RMSD) between two coordinate sets

Parameters:
  • coords1 (numpy.ndarray) – First coordinate set

  • coords2 (numpy.ndarray) – Second coordinate set

Returns:

RMSD value

Return type:

float

taco_format.center_of_mass(positions, masses=None)

Calculate the center of mass of a set of positions

Parameters:
  • positions (numpy.ndarray) – Atomic positions

  • masses (Optional[numpy.ndarray]) – Atomic masses (optional)

Returns:

Center of mass coordinates [x, y, z]

Return type:

numpy.ndarray

taco_format.cli_main()

CLI entry point for the installed ‘taco’ command

This function is called when the user runs ‘taco’ from the command line. It automatically gets the command line arguments from Python’s sys.argv.

taco_format.copy_frames(src_path, dst_path, start_frame, num_frames=None)

Copy frames from one TACO file to another

Parameters:
  • src_path (str) – Source file path

  • dst_path (str) – Destination file path

  • start_frame (int) – Index of the first frame to copy

  • num_frames (Optional[int]) – Number of frames to copy (None to copy all remaining frames)

taco_format.extract_atoms(src_path, dst_path, atom_indices)

Extract a subset of atoms from a TACO file

Parameters:
  • src_path (str) – Source file path

  • dst_path (str) – Destination file path

  • atom_indices (list[int]) – List of atom indices to extract

taco_format.extract_subset(positions, indices)

Extract a subset of atoms from positions

Parameters:
  • positions (numpy.ndarray) – Atomic positions

  • indices (list[int]) – Atom indices to extract

Returns:

Positions of the selected atoms

Return type:

numpy.ndarray

taco_format.get_file_info(path)

Get file information as a string

Parameters:

path (str) – Path to the TACO file

Returns:

Formatted string with file information

Return type:

str

taco_format.is_taco_file(path)

Check if a file is a TACO file

Parameters:

path (str) – Path to the file

Returns:

True if the file is a TACO file, False otherwise

Return type:

bool

taco_format.read(path, frame_index=0, read_as_trajectory=False, start_frame=0, end_frame=None)

Read ASE Atoms object(s) from a TACO file.

This function automatically determines whether the file contains a single structure or a trajectory and returns the appropriate object type.

Parameters:
  • path (str) – Path to the input file.

  • frame_index (int, optional) – Index of the frame to read for single structure (ignored for trajectory reading). Default is 0.

  • read_as_trajectory (bool, optional) – Force reading as a trajectory even for a single frame.

  • start_frame (int, optional) – Start index when reading a range of frames (inclusive). Default is 0.

  • end_frame (int, optional) – End index when reading a range of frames (exclusive). Default is None (read all frames).

Returns:

ASE Atoms object or list of ASE Atoms objects depending on the file content and parameters.

Return type:

Union[ase.Atoms, list[ase.Atoms]]

Example reading a single structure:

>>> import taco_format
>>> # This returns a single Atoms object if the file contains one frame
>>> # or a list of Atoms objects if it contains multiple frames
>>> atoms = taco_format.read('test.taco')

Example reading a specific frame from a trajectory:

>>> # Read the 10th frame (index 9) from a trajectory
>>> atoms = taco_format.read('traj.taco', frame_index=9)

Example forcing trajectory output:

>>> # Force reading as a trajectory even for a single frame
>>> traj = taco_format.read('single_frame.taco', read_as_trajectory=True)
>>> # traj will be a list with a single Atoms object

Example reading a range of frames:

>>> # Read frames 100 through 999
>>> frames = taco_format.read('traj.taco', start_frame=100, end_frame=1000)
taco_format.read_frame_range(path, start_frame, end_frame)

Read a range of frames from a TACO file.

Parameters:
  • path (str) – Path to the input file.

  • start_frame (int) – Index of the first frame to read.

  • end_frame – Index of the frame after the last one to read (exclusive).

  • end_frame – int

Returns:

List of ASE Atoms objects.

Return type:

list[ase.Atoms]

Example reading a range of frames:

>>> import taco_format
>>> # Read frames 10 through 19
>>> atoms_list = taco_format.read_frame_range('traj.taco', 10, 20)
taco_format.read_frames(path, frame_indices)

Read specific frames from a TACO file.

Parameters:
  • path (str) – Path to the input file

  • frame_indices (list[int]) – List of frame indices to read

Returns:

List of ASE Atoms objects

Return type:

list[ase.Atoms]

Example reading specific frames:

>>> import taco_format
>>> # Read frames 0, 10, 20, 30, 40
>>> frames = taco_format.read_frames('traj.taco', [0, 10, 20, 30, 40])
taco_format.run_cli(args)

Run the TACO CLI with the given arguments

This function provides access to the same CLI functionality as the standalone taco binary, but callable from Python.

Parameters:

args – List of command line arguments (same as you would pass to the taco binary) The first argument should be the command name (e.g., “info”, “extract”, etc.)

Example

>>> import taco_format
>>> taco_format.run_cli(["info", "trajectory.taco"])
>>> taco_format.run_cli(["extract", "input.taco", "output.taco", "--start", "0", "--end", "100"])
taco_format.write(path, atoms, time_step=0.001, full_frame_interval=10, compression_level=3, lossless=False)

Write ASE Atoms object(s) to a TACO file.

This function automatically detects whether it’s dealing with a single Atoms object or a list of Atoms objects and handles them accordingly.

Parameters:
  • path (str) – Path to the output file.

  • atoms (Union[ase.Atoms, list[ase.Atoms]]) – ASE Atoms object or list of ASE Atoms objects.

  • time_step (float, optional) – Time step between frames in picoseconds. Default is 0.001 (1 fs).

  • full_frame_interval (int, optional) – Interval for storing full frames (non-delta encoded) to prevent error accumulation. Default is 100.

  • compression_level (int, optional) – Compression level (1-22). Default is 3.

  • lossless (bool, optional) – Use lossless compression for all data. Default is False.

Example with a single Atoms object:

>>> import ase
>>> import taco_format
>>> atoms = ase.Atoms('H2O')
>>> taco_format.write('test.taco', atoms, compression_level=5)

Example with a trajectory:

>>> import ase
>>> import taco_format
>>> atoms_list = [ase.Atoms('H2O') for _ in range(10)]
>>> taco_format.write('traj.taco', atoms_list, time_step=0.002,
                     full_frame_interval=50,
                     compression_level=5)
taco_format.write_frames(path, atoms_list, time_step=0.001, full_frame_interval=10, compression_level=3, lossless=False)

Write a list of ASE Atoms objects to a TACO file.

Parameters:
  • path (str) – Path to the output file.

  • atoms_list (list[ase.Atoms]) – List of ASE Atoms objects.

  • time_step (float, optional) – Time step between frames in picoseconds. Default is 0.001 (1 fs).

  • full_frame_interval (int, optional) – Interval for storing full frames. Default is 100.

  • compression_level (int, optional) – Compression level (1-22). Default is 3.

  • lossless (bool, optional) – Use lossless compression for all data. Default is False.

Example:

>>> import taco_format
>>> atoms_list = [...] # List of ASE Atoms objects
>>> taco_format.write_frames('traj.taco', atoms_list)