Chemfiles.DataBufferType

Data (either binary or string) used for in-memory reading of trajectory

DataBuffer = Union{AbstractString,Vector{UInt8},Nothing}
Chemfiles.AtomType

An Atom is a particle in the current Frame.

An atom stores the following atomic properties:

  • atom name
  • atom type
  • atom mass
  • atom charge

The atom name is usually an unique identifier ("H1", "C_a") while the atom type will be shared among all particles of the same type: "H", "Ow", "CH3".

Chemfiles.AtomMethod
Atom(frame::Frame, index::Integer) -> Atom

Get a copy of the atom at the given index from a frame.

Chemfiles.AtomMethod
Atom(name::String) -> Atom

Create an atom with the given name and set the atom type to be the same as name.

Chemfiles.AtomMethod
Atom(topology::Topology, index::Integer) -> Atom

Get a copy of the atom at the given index from a topology.

Chemfiles.BondOrderType

Possible bond orders in Chemfiles:

  • Chemfiles.UnknownBond: when the bond order is not specified
  • Chemfiles.SingleBond: for single bonds
  • Chemfiles.DoubleBond: for double bonds
  • Chemfiles.TripleBond: for triple bonds
  • Chemfiles.QuadrupleBond: for quadruple bonds (present in some metals)
  • Chemfiles.QuintupletBond: for qintuplet bonds (present in some metals)
  • Chemfiles.AmideBond: for amide bonds
  • Chemfiles.AromaticBond: for aromatic bonds
Chemfiles.CellShapeType

The possible shape for an unit cell are:

  • Chemfiles.Orthorhombic for unit cells with the three angles are 90°.
  • Chemfiles.Triclinic for unit cells where the three angles may not be 90°.
  • Chemfiles.Infinite for unit cells without boundaries.
Chemfiles.FormatMetadataType

Metadata associated with one of the format Chemfiles can read/write

  • name::String: Name of the format

  • extension::Union{Nothing, String}: Extension associated with the format, or nothing if there is no extension associated.

  • description::String: Extended user-facing description of the format

  • reference::String: URL pointing to the format definition/reference

  • read::Bool: Is reading files in this format implemented?

  • write::Bool: Is writing files in this format implemented?

  • memory::Bool: Does this format support in-memory IO?

  • positions::Bool: Does this format support storing atomic positions?

  • velocities::Bool: Does this format support storing atomic velocities?

  • unit_cell::Bool: Does this format support storing unit cell information?

  • atoms::Bool: Does this format support storing atom names or types?

  • bonds::Bool: Does this format support storing bonds between atoms?

  • residues::Bool: Does this format support storing residues?

Chemfiles.FrameType

A Frame holds data for one step of a simulation. As not all formats provide all the types of information, some fields may be initialized to a default value. A Frame may contain the following data:

  • Positions for all the atoms in the system;
  • Velocities for all the atoms in the system;
  • The Topology of the system;
  • The UnitCell of the system.
Chemfiles.ResidueType
Residue(name::String) -> Residue
Residue(name::String, id) -> Residue

Create a new residue with the given name and optional residue identifier id.

Chemfiles.ResidueType

A Residue is a group of atoms belonging to the same logical unit. They can be small molecules, amino-acids in a protein, monomers in polymers, etc.

Chemfiles.ResidueMethod
Residue(topology::Topology, index::Integer) -> Residue

Get a copy of the residue at index from a topology.

The residue index in the topology is not always the same as the residue identifier.

Chemfiles.SelectionType

A Selection is used to select a group of atoms. Examples of selections are "name H" and "(x < 45 and name O) or name C". See the full documentation <http://chemfiles.org/chemfiles/latest/selections.html>_ for more information about the selection language.

Chemfiles.SelectionMethod
Selection(selection::AbstractString) -> Selection

Create a Selection from a selection string.

Chemfiles.TopologyType

A Topology describes the organisation of the particles in the system: what their names are, how they are bonded together, etc. A Topology is a list of Atoms in the system, together with the list of bonds between the atoms.

Chemfiles.TopologyMethod
Topology(frame::Frame) -> Topology

Get a copy of the Topology of the given frame.

Chemfiles.TrajectoryType
Trajectory(path::AbstractString) -> Trajectory
Trajectory(path::AbstractString, mode::Char) -> Trajectory
Trajectory(
    path::AbstractString,
    mode::Char,
    format::AbstractString
) -> Trajectory

Opens the trajectory file at the given path.

The opening mode can be 'r' for read, 'w' for write or 'a' for append, and defaults to 'r'. The optional format parameter give a specific file format to use when opening the file.

Chemfiles.TrajectoryType

A Trajectory represents a simulation file on the hard drive. It can read or write one or many Frames to this file. The file format can be automatically determined from the extention, or manually specified. Writing to a Trajectory is buffered, which means that one needs to close() the trajectory and flush the buffer before being able to read the file again.

Chemfiles.TrajectoryMethod
Trajectory(f::Function, args...) -> Any

Apply the function f to the result of Trajectory(args...) and close the resulting trajectory upon completion, similar to open(f, args...).

Chemfiles.UnitCellType
UnitCell(lengths::Vector{Float64}) -> UnitCell
UnitCell(
    lengths::Vector{Float64},
    angles::Vector{Float64}
) -> UnitCell

Create an UnitCell from the given lengths and angles.

Chemfiles.UnitCellType

A UnitCell describes the bounding box of a system. It is represented by a 3x3 matrix containing the base vectors a, b and c.

Chemfiles.UnitCellMethod
UnitCell(frame::Frame) -> UnitCell

Get a copy of the UnitCell of a frame.

Chemfiles.UnitCellMethod
UnitCell(matrix::Matrix{Float64}) -> UnitCell

Create an UnitCell from the given 3x3 cell matrix.

AtomsBase.atomic_numberMethod
atomic_number(atom::Atom) -> UInt64

Get the atomic number of an atom from the atom type.

If the atomic number can not be found, returns 0.

Base.angleMethod
angle(
    frame::Frame,
    i::Integer,
    j::Integer,
    k::Integer
) -> Float64

Calculate the angle made by three atoms.

Base.closeMethod
close(trajectory::Trajectory)

Close a trajectory. This function flushes any buffer content to the hard drive, and frees the associated memory. Necessary when running on the REPL to finish writing.

Base.containsMethod
contains(residue::Residue, index::Integer) -> Bool

Check if the atom at the given index is in the residue.

Base.convertMethod
convert(
    _::Type{AtomsBase.AbstractSystem},
    frame::Frame
) -> AtomsBase.FlexibleSystem

Convert a Chemfiles Frame to an AtomsBase AbstractSystem

Base.convertMethod
convert(
    _::Type{AtomsBase.FlexibleSystem},
    frame::Frame
) -> AtomsBase.FlexibleSystem

Convert a Chemfiles Frame to an AtomsBase FlexibleSystem

Base.deepcopyMethod
deepcopy(atom::Atom) -> Atom

Make a deep copy of an atom.

Base.deepcopyMethod
deepcopy(frame::Frame) -> Frame

Make a deep copy of a Frame.

Base.deepcopyMethod
deepcopy(residue::Residue) -> Residue

Make a deep copy of a residue.

Base.deepcopyMethod
deepcopy(selection::Selection) -> Selection

Make a deep copy of a selection.

Base.deepcopyMethod
deepcopy(topology::Topology) -> Topology

Make a deep copy of a topology.

Base.deepcopyMethod
deepcopy(cell::UnitCell) -> UnitCell

Make a deep copy of a cell.

Base.fullnameMethod
fullname(atom::Atom) -> String

Get the full name of an atom from the atom type.

For example, the full name of an atom with type "He" is "Helium".

Base.isopenMethod
isopen(trajectory::Trajectory) -> Bool

Check if the trajectory is open.

Base.lengthMethod
length(frame::Frame) -> Int64

Get the number of atoms in the frame.

Base.lengthMethod
length(trajectory::Trajectory) -> UInt64

Get the number of steps (the number of frames) in a trajectory.

Base.read!Method
read!(trajectory::Trajectory, frame::Frame)

Read the next step of the trajectory data into the preexisting Frame structure.

Base.readMethod
read(trajectory::Trajectory) -> Frame

Read the next step of the trajectory, and return the corresponding Frame.

Base.resize!Method
resize!(frame::Frame, natoms::Integer)

Resize the positions and the velocities in the frame, to make space for natoms atoms. This function may invalidate any pointer to the positions or the velocities if the new size is bigger than the old one. In all the cases, previous data is conserved. This function conserve the presence or absence of velocities.

Base.resize!Method
resize!(topology::Topology, size::Integer)

Resize the topology to hold natoms atoms. If the new number of atoms is bigger than the current number, new atoms will be created with an empty name and type.

If it is lower than the current number of atoms, the last atoms will be removed, together with the associated bonds, angles and dihedrals.

Base.sizeMethod
size(frame::Frame) -> Int64

Get the number of atoms in the frame.

Base.sizeMethod
size(residue::Residue) -> Int64

Get the number of atoms in a residue.

Base.sizeMethod
size(selection::Selection) -> UInt64

Get the size of a selection, i.e. the number of atoms we are selecting together.

Base.sizeMethod
size(topology::Topology) -> Int64

Get the Topology size, i.e. the current number of atoms.

Base.sizeMethod
size(trajectory::Trajectory) -> UInt64

Get the number of steps (the number of frames) in a trajectory.

Base.stepMethod
step(frame::Frame) -> UInt64

Get the frame step, i.e. the frame number in the trajectory.

Base.take!Method
take!(trajectory::Trajectory) -> Vector{UInt8}

Obtain the memory buffer of a in-memory trajectory writer as a Vector{UInt8}, without copying.

If more data is written to the trajectory, the buffer is invalidated.

Base.viewMethod
view(frame::Frame, index::Integer) -> Atom

Get the Atom at the given index of the frame without creating a copy.

Warning

This function can lead to unefined behavior when keeping the returned Atom around. Whith code like this:

frame = Frame()
resize!(frame, 3)

atom = @view frame[0]
resize!(frame, 4)

atom contains a pointer to memory owned by frame, but this pointer has been invalidated when resizing the frame. Using atom after the second call to resize! might trigger undefined behavior (segmentation fault in the best case, silent data corruption in the worst case).

Base.viewMethod
view(topology::Topology, index::Integer) -> Atom

Get the Atom at the given index of the topology without creating a copy.

Warning

This function can lead to unefined behavior when keeping the returned Atom around. See Base.view(frame::Frame, index::Integer) for more information on this issue.

Base.writeMethod
write(trajectory::Trajectory, frame::Frame)

Write the given frame to the trajectory.

Chemfiles.MemoryTrajectoryFunction
MemoryTrajectory(format::AbstractString) -> Trajectory
MemoryTrajectory(
    format::AbstractString,
    data::Union{Nothing, AbstractString, Vector{UInt8}}
) -> Trajectory

Open a trajectory that read in-memory data as if it were a formatted file.

The format of the data is mandatory, and must match one of the known formats. If data is nothing, this function creates a memory writer, the resulting data can be retrieved with buffer. Else, this function creates a memory reader using the given data.

Chemfiles.add_atom!Function
add_atom!(
    frame::Frame,
    atom::Atom,
    position::Vector{Float64}
)
add_atom!(
    frame::Frame,
    atom::Atom,
    position::Vector{Float64},
    velocity::Vector{Float64}
)

Add an atom and the corresponding position and velocity data to a frame.

Chemfiles.add_atom!Method
add_atom!(residue::Residue, index::Integer)

Add the atom at the given index in the residue.

Chemfiles.add_atom!Method
add_atom!(topology::Topology, atom::Atom)

Add an atom at the end of a topology.

Chemfiles.add_bond!Function
add_bond!(topology::Topology, i::Integer, j::Integer)
add_bond!(topology::Topology, i::Integer, j::Integer, order)

Add a bond between the atoms i and j in the topology, optionaly setting the bond order.

Chemfiles.add_bond!Function
add_bond!(frame::Frame, i::Integer, j::Integer)
add_bond!(frame::Frame, i::Integer, j::Integer, order)

Add an additional bond to the Frame's Topology.

Chemfiles.add_configurationMethod
add_configuration(path::String)

Read configuration data from the file at path.

By default, chemfiles reads configuration from any file name .chemfilesrc in the current directory or any parent directory. This function can be used to add data from another configuration file.

This function will fail if there is no file at path, or if the file is incorectly formatted. Data from the new configuration file will overwrite any existing data.

Chemfiles.add_residue!Method
add_residue!(frame::Frame, residue::Residue)

Add a residue to the Frame's Topology.

Chemfiles.add_residue!Method
add_residue!(topology::Topology, residue::Residue)

Add a copy of residue to this topology.

The residue id must not already be in the topology, and the residue must contain only atoms that are not already in another residue.

Chemfiles.add_velocities!Method
add_velocities!(frame::Frame)

Add velocities to this frame. The storage is initialized with the result of size(frame) as the number of atoms. If the frame already has velocities, this does nothing.

Chemfiles.anglesMethod
angles(topology::Topology) -> Matrix{UInt64}

Get the angles in the topology, in a 3 x angles_count(topology) array.

Chemfiles.anglesMethod
angles(cell::UnitCell) -> Vector{Float64}

Get the three cell angles in degrees.

Chemfiles.angles_countMethod
angles_count(topology::Topology) -> Int64

Get the number of angles in the topology.

Chemfiles.are_linkedMethod
are_linked(
    topology::Topology,
    first::Residue,
    second::Residue
) -> Bool

Check if the two residues first and second from the topology are linked together, i.e. if there is a bond between one atom in the first residue and one atom in the second one.

Chemfiles.atomsMethod
atoms(residue::Residue) -> Vector{UInt64}

Get the atoms in a $residue$. This function returns a list of indexes.

Chemfiles.bond_orderMethod
bond_order(
    topology::Topology,
    i::Integer,
    j::Integer
) -> Chemfiles.BondOrder

Get the BondOrder for the bond between atoms i and j in the topology.

Chemfiles.bond_ordersMethod
bond_orders(
    topology::Topology
) -> Vector{Chemfiles.BondOrder}

Get the BondOrder for all the bonds in the topology.

Chemfiles.bondsMethod
bonds(topology::Topology) -> Matrix{UInt64}

Get the bonds in the topology, in a 2 x bonds_count(topology) array.

Chemfiles.bonds_countMethod
bonds_count(topology::Topology) -> Int64

Get the number of bonds in the topology.

Chemfiles.chargeMethod
charge(atom::Atom) -> Float64

Get the charge of an atom in number of the electron charge e.

Chemfiles.clear_bonds!Method
clear_bonds!(frame::Frame)

Remove all bonds, angles and dihedral angles from the Frame's Topology.

Chemfiles.clear_bonds!Method
clear_bonds!(topology::Topology)

Remove all bonds, angles and dihedral angles from this Topology.

Chemfiles.count_residuesMethod
count_residues(topology::Topology) -> Int64

Get the number of residues in the topology.

Chemfiles.covalent_radiusMethod
covalent_radius(atom::Atom) -> Float64

Get the covalent radius of an atom from the atom type.

If the radius can not be found, returns 0.

Chemfiles.dihedralMethod
dihedral(
    frame::Frame,
    i::Integer,
    j::Integer,
    k::Integer,
    m::Integer
) -> Float64

Calculate the dihedral (torsional) angle made by four unbranched atoms.

Chemfiles.dihedralsMethod
dihedrals(topology::Topology) -> Matrix{UInt64}

Get the dihedral angles in the topology, in a 4 x dihedrals_count(topology) array.

Chemfiles.dihedrals_countMethod
dihedrals_count(topology::Topology) -> Int64

Get the number of dihedral angles in the topology.

Chemfiles.distanceMethod
distance(frame::Frame, i::Integer, j::Integer) -> Float64

Calculate the distance between two atoms.

Chemfiles.evaluateMethod
evaluate(selection::Selection, frame::Frame) -> Vector{Any}

Evaluate a selection on a given frame. This function return a list of indexes or tuples of indexes of atoms in the frame matching the selection.

Chemfiles.format_listMethod
format_list() -> Vector{Chemfiles.FormatMetadata}

Get the full list of formats supported by Chemfiles, and associated metadata

Chemfiles.guess_bonds!Method
guess_bonds!(frame::Frame)

Guess the bonds, angles, and dihedrals in the frame using a distance criteria.

Chemfiles.guess_formatMethod
guess_format(path::String) -> String

Get the format that chemfiles would use to read a file at the given path.

The format is mostly guessed from the path extension, chemfiles only tries to read the file to distinguish between CIF and mmCIF files. Opening the file using the returned format string might still fail. For example, it will fail if the file is not actually formatted according to the guessed format; or the format/compression combination is not supported (e.g. XTC / GZ will not work since the XTC reader does not support compressed files).

The returned format is represented in a way compatible with the various Trajectory constructors, i.e. "<format name> [/ <compression>]", where compression is optional.

Chemfiles.idMethod
id(residue::Residue) -> Int64

Get the identifier of a residue in the initial topology.

Chemfiles.impropersMethod
impropers(topology::Topology) -> Matrix{UInt64}

Get the improper angles in the topology, in a 4 x impropers_count(topology) array.

Chemfiles.impropers_countMethod
impropers_count(topology::Topology) -> Int64

Get the number of improper angles in the topology.

Chemfiles.last_errorMethod
last_error() -> String

Get the last error message from the chemfiles runtime.

Chemfiles.lengthsMethod
lengths(cell::UnitCell) -> Vector{Float64}

Get the three cell lengths in angstroms.

Chemfiles.list_propertiesMethod
list_properties(atom::Atom) -> Vector{String}

Get the names of all properties associated with an atom.

Chemfiles.list_propertiesMethod
list_properties(frame::Frame) -> Vector{String}

Get the names of all properties associated with a frame.

Chemfiles.list_propertiesMethod
list_properties(residue::Residue) -> Vector{String}

Get the names of all properties associated with a residue.

Chemfiles.massMethod
mass(atom::Atom) -> Float64

Get the mass of an atom in atomic mass units.

Chemfiles.matrixMethod
matrix(cell::UnitCell) -> Matrix{Float64}

Get the cell matricial representation, i.e. the representation of the three base vectors as:

    | a_x   b_x   c_x |
    |  0    b_y   c_y |
    |  0     0    c_z |
Chemfiles.nameMethod
name(atom::Atom) -> String

Get the name of an atom.

Chemfiles.nameMethod
name(residue::Residue) -> String

Get the name of a residue.

Chemfiles.out_of_planeMethod
out_of_plane(
    frame::Frame,
    i::Integer,
    j::Integer,
    k::Integer,
    m::Integer
) -> Float64

Calculate the out-of-plane (improper) angle made by four atoms.

Chemfiles.pathMethod
path(trajectory::Trajectory) -> String

Get the path used to open a trajectory.

Chemfiles.positionsMethod
positions(frame::Frame) -> Chemfiles.ChemfilesArray

Get the positions in a Frame as an array. The positions are readable and writable from this array. If the frame is resized (by writing to it, or calling resize!), the array is invalidated.

Chemfiles.properties_countMethod
properties_count(frame::Frame) -> Int64

Get the number of properties associated with a frame.

Chemfiles.properties_countMethod
properties_count(residue::Residue) -> Int64

Get the number of properties associated with a residue.

Chemfiles.propertyMethod
property(atom::Atom, name::String) -> Any

Get a named property for the given atom.

Chemfiles.propertyMethod
property(frame::Frame, name::String) -> Any

Get a named property for the given atom.

Chemfiles.propertyMethod
property(residue::Residue, name::String) -> Any

Get a named property for the given residue.

Chemfiles.read_step!Method
read_step!(
    trajectory::Trajectory,
    step::Integer,
    frame::Frame
)

Read the given step of the trajectory into an preexisting Frame structure.

Chemfiles.read_stepMethod
read_step(trajectory::Trajectory, step::Integer) -> Frame

Read the given step of the trajectory, and return the corresponding Frame.

Chemfiles.remove_atom!Method
remove_atom!(frame::Frame, index::Integer)

Remove the atom at index from the frame.

This function modifies all the atoms indexes after index, and invalidates any array obtained using positions or velocities.

Chemfiles.remove_atom!Method
remove_atom!(topology::Topology, index::Integer)

Remove the atom at the given index from a topology.

Chemfiles.remove_bond!Method
remove_bond!(frame::Frame, i::Integer, j::Integer)

Remove a bond from the Frame's Topology.

Chemfiles.remove_bond!Method
remove_bond!(topology::Topology, i::Integer, j::Integer)

Remove any existing bond between the atoms i and j in the topology.

Chemfiles.residue_for_atomMethod
residue_for_atom(
    topology::Topology,
    index::Integer
) -> Union{Nothing, Residue}

Get a copy of the residue containing the atom at index in the topology.

This function will return nothing if the atom is not in a residue, or if the index is bigger than the number of atoms in the topology.

Chemfiles.selection_stringMethod
selection_string(selection::Selection) -> String

Get the selection string used to create a given selection.

Chemfiles.set_angles!Method
set_angles!(cell::UnitCell, angles::Vector{Float64})

Set the cell angles to the given values. The angles should be in degrees.

Chemfiles.set_cell!Method
set_cell!(frame::Frame, cell::UnitCell)

Set the cell associated with a frame.

Chemfiles.set_cell!Method
set_cell!(trajectory::Trajectory, cell::UnitCell)

Set the cell associated with a trajectory. This cell will be used when reading and writing the file, replacing any unit cell in the file.

Chemfiles.set_charge!Method
set_charge!(atom::Atom, charge)

Set the charge of an atom to charge.

The charge must be in number of the electron charge e.

Chemfiles.set_lengths!Method
set_lengths!(cell::UnitCell, lengths::Vector{Float64})

Set the cell lengths to the given values. The lengths should be in angstroms.

Chemfiles.set_mass!Method
set_mass!(atom::Atom, mass)

Set the mass of an atom to mass.

The mass must be in atomic mass units.

Chemfiles.set_name!Method
set_name!(atom::Atom, name::String)

Set the name of an atom to name.

Chemfiles.set_property!Method
set_property!(
    atom::Atom,
    name::String,
    value::Union{Bool, Float64, String, Vector{Float64}}
)

Set a named property for the given atom.

Chemfiles.set_property!Method
set_property!(frame::Frame, name::String, value)

Set a named property for the given Frame.

Chemfiles.set_property!Method
set_property!(
    residue::Residue,
    name::String,
    value::Union{Bool, Float64, String, Vector{Float64}}
)

Set a named property for the given residue.

Chemfiles.set_shape!Method
set_shape!(cell::UnitCell, shape::CellShape)

Set the cell shape to the given shape.

Chemfiles.set_step!Method
set_step!(frame::Frame, step::Integer)

Set the frame step to step.

Chemfiles.set_topology!Function
set_topology!(trajectory::Trajectory, path::AbstractString)
set_topology!(
    trajectory::Trajectory,
    path::AbstractString,
    format::AbstractString
)

Set the Topology associated with a trajectory by reading the first frame of the file at path; and extracting the topology of this frame. The optional format parameter can be used to specify the file format.

Chemfiles.set_topology!Method
set_topology!(frame::Frame, topology::Topology)

Set the topology associated with a frame.

Chemfiles.set_topology!Method
set_topology!(trajectory::Trajectory, topology::Topology)

Set the Topology associated with a trajectory. This topology will be used when reading and writing the file, replacing any topology in the file.

Chemfiles.set_type!Method
set_type!(atom::Atom, type::String)

Set the type of an atom to type.

Chemfiles.set_warning_callbackMethod
set_warning_callback(callback::Function)

Set the global warning callback to be used for each warning event.

The callback function must take a String and return nothing.

Chemfiles.shapeMethod
shape(cell::UnitCell) -> CellShape

Get the cell shape, as a CellShape value.

Chemfiles.typeMethod
type(atom::Atom) -> String

Get the type of an atom.

Chemfiles.vdw_radiusMethod
vdw_radius(atom::Atom) -> Float64

Get the van der Waals radius of an atom from the atom type.

If the radius can not be found, this function returns 0.

Chemfiles.velocitiesMethod
velocities(frame::Frame) -> Chemfiles.ChemfilesArray

Get the velocities in a Frame as an array. The velocities are readable and writable from this array. If the frame is resized (by writing to it, or calling resize!), the array is invalidated.

If the frame do not have velocity, this function will error. You can use add_velocities! to add velocities to a frame before calling this function.

Chemfiles.volumeMethod
volume(cell::UnitCell) -> Float64

Get the unit cell volume.

Chemfiles.wrap!Method
wrap!(
    cell::UnitCell,
    vector::Vector{Float64}
) -> Vector{Float64}

Wrap a vector in the unit cell.