Biomolecular systems

Systems

BiochemicalAlgorithms.SystemType
mutable struct System{T} <: AbstractAtomContainer{T}

Mutable representation of a biomolecular system.

Fields

  • name::String
  • properties::Properties
  • flags::Flags

Constructors

System(name::AbstractString = "", properties::Properties = Properties(), flags::Flags = Flags())

Creates a new and empty System{Float32}.

System{T}(name::AbstractString = "", properties::Properties = Properties(), flags::Flags = Flags())

Creates a new and empty System{T}.

source
BiochemicalAlgorithms.parent_systemFunction
parent_system(::Atom)
parent_system(::Bond)
parent_system(::Chain)
parent_system(::SecondaryStructure)
parent_system(::Fragment)
parent_system(::Molecule)
parent_system(::System)

Returns the System{T} containing the given object. Alias for Base.parent.

source
Base.parentMethod
parent(::Atom)
parent(::Bond)
parent(::Chain)
parent(::SecondaryStructure)
parent(::Fragment)
parent(::Molecule)
parent(::System)

Returns the System{T} containing the given object.

source

System components

BiochemicalAlgorithms.SystemComponentTableType
struct SystemComponentTable{T, C <: AbstractSystemComponent{T}} <: AbstractSystemComponentTable{T}

Tables.jl-compatible system component table for a specific System{T} and system component C (e.g., Atom{T}, Bond{T}, etc.).

source
BiochemicalAlgorithms.get_propertyFunction
get_property(
    ac::AbstractSystemComponent,
    key::Symbol
) -> Any

Returns the property associated with the given key in ac.

source
get_property(
    ac::AbstractSystemComponent,
    key::Symbol,
    default
) -> Any

Returns the property associated with the given key in ac. If no such property exists, returns default.

source
BiochemicalAlgorithms.has_flagFunction
has_flag(ac::AbstractSystemComponent, flag::Symbol) -> Any

Returns a Bool indicating whether the given system component has the given flag.

source
BiochemicalAlgorithms.has_propertyFunction
has_property(
    ac::AbstractSystemComponent,
    key::Symbol
) -> Any

Returns a Bool indicating whether the given system component has the given property.

source
BiochemicalAlgorithms.revalidate_indices!Function
revalidate_indices!(::SystemComponentTable)
revalidate_indices!(::SystemComponentTableCol)

Removes remnants of previously removed system components from the given table or table column.

source
Base.sortMethod
sort(::SystemComponentTable)

Returns a copy of the given table, sorted by idx (default) or according to the given keyword arguments.

Supported keyword arguments

Same as Base.sort

source
Base.sort!Method
sort!(::SystemComponentTable)

In-place variant of sort.

Note

Only the given table is modified, not the underlying system!

source

Atoms

BiochemicalAlgorithms.AtomType
Atom{T} <: AbstractSystemComponent{T}

Mutable representation of an individual atom in a system.

Public fields

  • idx::Int
  • number::Int
  • element::ElementType
  • name::String
  • atom_type::String
  • r::Vector3{T}
  • v::Vector3{T}
  • F::Vector3{T}
  • formal_charge::Int
  • charge::T
  • radius::T

Private fields

  • properties::Properties
  • flags::Flags
  • frame_id::Int
  • molecule_idx::MaybeInt
  • chain_idx::MaybeInt
  • fragment_idx::MaybeInt

Constructors

Atom(
    ac::AbstractAtomContainer{T}
    number::Int,
    element::ElementType;
    # keyword arguments
    name::AbstractString = "",
    atom_type::AbstractString = "",
    r::Vector3{T} = Vector3{T}(0, 0, 0),
    v::Vector3{T} = Vector3{T}(0, 0, 0),
    F::Vector3{T} = Vector3{T}(0, 0, 0),
    formal_charge::Int = 0,
    charge::T = zero(T),
    radius::T = zero(T),
    properties::Properties = Properties(),
    flags::Flags = Flags(),
    frame_id::Int = 1
    molecule_idx::MaybeInt = nothing,
    chain_idx::MaybeInt = nothing,
    fragment_idx::MaybeInt = nothing
)

Creates a new Atom{T} in the given atom container.

Atom(
    number::Int,
    element::ElementType;
    kwargs...
)

Creates a new Atom{Float32} in the default system. Supports the same keyword arguments as above.

source
BiochemicalAlgorithms.AtomTableType
AtomTable{T} <: AbstractSystemComponentTable{T}

Tables.jl-compatible representation of system atoms (or a subset thereof). Atom tables can be generated using atoms or filtered from other atom tables (via Base.filter).

Public columns

  • idx::AbstractVector{Int}
  • number::AbstractVector{Int}
  • element::AbstractVector{ElementType}
  • name::AbstractVector{String}
  • atom_type::AbstractVector{String}
  • r::AbstractVector{Vector3{T}}
  • v::AbstractVector{Vector3{T}}
  • F::AbstractVector{Vector3{T}}
  • formal_charge::AbstractVector{Int}
  • charge::AbstractVector{T}
  • radius::AbstractVector{T}

Private columns

  • properties::AbstractVector{Properties}
  • flags::AbstractVector{Flags}
  • frame_id::AbstractVector{Int}
  • molecule_idx::AbstractVector{MaybeInt}
  • chain_idx::AbstractVector{MaybeInt}
  • fragment_idx::AbstractVector{MaybeInt}
source
BiochemicalAlgorithms.atom_by_idxFunction
atom_by_idx(
    sys::System{T} = default_system(),
    idx::Int
) -> Atom{T}

Returns the Atom{T} associated with the given idx in sys. Throws a KeyError if no such atom exists.

source
BiochemicalAlgorithms.atom_by_nameFunction
atom_by_name(
    ac::AbstractAtomContainer{T} = default_system(),
    name::AbstractString
) -> Union{Nothing, Atom{T}}

Returns the first Atom{T} associated with the given name in ac. Returns nothing if no such atom exists.

Supported keyword arguments

  • frame_id::MaybeInt = 1: Any value other than nothing limits the result to atoms matching this frame ID.
source
BiochemicalAlgorithms.atomsFunction
atoms(::Chain)
atoms(::Fragment)
atoms(::Molecule)
atoms(::System = default_system())

Returns an AtomTable{T} containing all atoms of the given atom container.

Supported keyword arguments

  • frame_id::MaybeInt = 1
  • molecule_idx::Union{MaybeInt, Some{Nothing}} = nothing
  • chain_idx::Union{MaybeInt, Some{Nothing}} = nothing
  • fragment_idx::Union{MaybeInt, Some{Nothing}} = nothing

All keyword arguments limit the results to atoms matching the given IDs. Keyword arguments set to nothing are ignored. You can use Some(nothing) to explicitly filter for ID values of nothing.

source
atoms(::ChainTable)
atoms(::FragmentTable)
atoms(::MoleculeTable)

Returns an AtomTable{T} containing all atoms of the given table.

source
atoms(
    substruct::Substructure{T, A} where A<:AbstractAtomContainer{T};
    frame_id,
    molecule_idx,
    chain_idx,
    secondary_structure_idx,
    fragment_idx
) -> SystemComponentTable{T, C} where {T, C<:Atom{T}}

Returns an AtomTable for all of the given system's atoms matching the given criteria (value or missing). Fields given as nothing are ignored. The returned table contains all public and private atom fields.

source
atoms(
    ab::AbstractAtomBijection
) -> Tuple{AtomTable, AtomTable}

Returns the tuple of atom tables represented by this bijection.

source
BiochemicalAlgorithms.is_bound_toFunction
is_bound_to(a1::Atom, a2::Atom) -> Bool

Decides if two atoms are bound to each other. Hydrogen bonds (has_flag(bond, :TYPE__HYDROGEN)) are ignored.

source
BiochemicalAlgorithms.is_geminalFunction
is_geminal(a1::Atom, a2::Atom) -> Union{Missing, Bool}

Decides if two atoms are geminal.

Two atoms are geminal if they do not share a common bond but both have a bond to a third atom. For example the two hydrogen atoms in water are geminal. Hydrogen bonds (has_flag(bond, :TYPE__HYDROGEN)) are ignored.

source
BiochemicalAlgorithms.is_vicinalFunction
is_vicinal(a1::Atom, a2::Atom) -> Bool

Decides if two atoms are vicinal.

Two atoms are vicinal if they are separated by three bonds (1-4 position). Hydrogen bonds (has_flag(bond, :TYPE__HYDROGEN)) are ignored.

source
BiochemicalAlgorithms.natomsFunction
natoms(::Chain)
natoms(::Fragment)
natoms(::Molecule)
natoms(::System = default_system())

Returns the number of atoms in the given atom container.

Supported keyword arguments

See atoms

source
natoms(::AtomTable)
natoms(::ChainTable)
natoms(::FragmentTable)
natoms(::MoleculeTable)

Returns the number of atoms in the given table.

source
BiochemicalAlgorithms.sort_atoms!Function
sort_atoms!(::System)

Sorts the atoms in the given system by idx (default) or according to the given keyword arguments.

Supported keyword arguments

Same as Base.sort

source
Base.delete!Method
delete!(::Atom)
delete!(::AtomTable)
delete!(::AtomTable, idx::Int)

Removes the given atom(s) and all associated bonds from the associated system.

source
Base.push!Method
push!(::Fragment{T},   ::Atom{T})
push!(::Molecule{T},   ::Atom{T})
push!(::System{T},     ::Atom{T})

Creates a copy of the given atom in the given atom container. The new atom is automatically assigned a new idx.

Supported keyword arguments

See atoms

source

Bonds

BiochemicalAlgorithms.BondType
Bond{T} <: AbstractAtomContainer{T}

Mutable representation of an individual bond in a system.

Public fields

  • idx::Int
  • a1::Int
  • a2::Int
  • order::BondOrderType

Private fields

  • properties::Properties
  • flags::Flags

Constructors

Bond(
    ac::AbstractAtomContainer{T} = default_system(),
    a1::Int,
    a2::Int,
    order::BondOrderType;
    # keyword arguments
    properties::Properties = Properties(),
    flags::Flags = Flags()
)

Creates a new Bond{T} in the given (atom container's) system.

Bond(
    a1::Atom{T},
    a2::Atom{T},
    order::BondOrderType;
    # keyword arguments
    properties::Properties = Properties(),
    flags::Flags = Flags()
)

Creates a new Bond{T} for the given atoms. Both atoms must belong to the same system.

source
BiochemicalAlgorithms.BondTableType
BondTable{T} <: AbstractSystemComponentTable{T}

Tables.jl-compatible representation of system bonds (or a subset thereof). Bond tables can be generated using bonds or filtered from other bond tables (via Base.filter).

Public columns

  • idx::AbstractVector{Int}
  • a1::AbstractVector{Int}
  • a2::AbstractVector{Int}
  • order::AbstractVector{BondOrderType}

Private columns

  • properties::AbstractVector{Properties}
  • flags::AbstractVector{Flags}
source
BiochemicalAlgorithms.bond_by_idxFunction
bond_by_idx(
    sys::System{T} = default_system(),
    idx::Int
) -> Bond{T}

Returns the Bond{T} associated with the given idx in sys. Throws a KeyError if no such bond exists.

source
BiochemicalAlgorithms.bondsFunction
bonds(::Atom)
bonds(::AtomTable)

Returns a BondTable{T} containing all bonds of the given atom(s).

source
bonds(::Chain)
bonds(::Fragment)
bonds(::Molecule)
bonds(::System = default_system())

Returns a BondTable{T} containing all bonds of the given atom container where at least one associated atom is contained in the same container.

Supported keyword arguments

See atoms

source
bonds(::ChainTable)
bonds(::FragmentTable)
bonds(::MoleculeTable)

Returns a BondTable{T} containing all bonds where at least one associated atom belongs to the given table.

source
BiochemicalAlgorithms.nbondsFunction
nbonds(::Atom)
nbonds(::AtomTable)

Returns the number of bonds of the given atom(s).

source
nbonds(::Chain)
nbonds(::Fragment)
nbonds(::Molecule)
nbonds(::System = default_system())

Returns the number of bonds in the given atom container where at least one associated atom is contained in the same container.

Supported keyword arguments

See atoms

source
nbonds(::BondTable)
nbonds(::ChainTable)
nbonds(::FragmentTable)
nbonds(::MoleculeTable)

Returns the number of bonds where at least one associated atom belongs to the given table.

source
BiochemicalAlgorithms.sort_bonds!Function
sort_bonds!(::System)

Sorts the bonds in the given system by idx (default) or according to the given keyword arguments.

Supported keyword arguments

Same as Base.sort

source
Base.delete!Method
delete!(::Bond)
delete!(::BondTable)
delete!(::BondTable, idx::Int)

Removes the given bond(s) from the associated system.

source
Base.push!Method
push!(::AbstractAtomContainer, ::Bond{T})

Creates a copy of the given bond in the system associated with the given atom container. The new bond is automatically assigned a new idx.

source

Molecules

BiochemicalAlgorithms.MoleculeVariantModule
@enumx MoleculeVariant begin
    None = 1
    Protein = 2
end

Enum representing variants of molecules

Example

julia> prot = Protein(System())
Molecule{Float32}: (idx = 1, name = "")

julia> isprotein(prot)
true

julia> prot.variant == MoleculeVariant.Protein
true
source

Molecules (all variants)

BiochemicalAlgorithms.MoleculeType
Molecule{T} <: AbstractAtomContainer{T}

Mutable representation of an individual molecule in a system.

Public fields

  • idx::Int
  • name::String

Private fields

  • variant::MoleculeVariantType
  • properties::Properties
  • flags::Flags

Constructors

Molecule(
    sys::System{T};
    # keyword arguments
    name::AbstractString = "",
    variant::MoleculeVariantType = MoleculeVariant.None,
    properties::Properties = Properties(),
    flags::Flags = Flags()
)

Creates a new Molecule{T} in the given system.

Molecule(;
    #keyword arguments
    name::AbstractString = "",
    variant::MoleculeVariantType = MoleculeVariant.None,
    properties::Properties = Properties(),
    flags::Flags = Flags()
)

Creates a new Molecule{Float32} in the default system.

source
BiochemicalAlgorithms.MoleculeTableType
MoleculeTable{T} <: AbstractSystemComponentTable{T}

Tables.jl-compatible representation of system molecules (or a subset thereof). Molecule tables can be generated using molecules or filtered from other molecule tables (via Base.filter).

Public columns

  • idx::AbstractVector{Int}
  • name::AbstractVector{String}

Private columns

  • variant::AbstractVector{MoleculeVariantType}
  • properties::AbstractVector{Properties}
  • flags::AbstractVector{Flags}
source
BiochemicalAlgorithms.molecule_by_idxFunction
molecule_by_idx(
    sys::System{T} = default_system(),
    idx::Int
) -> Molecule{T}

Returns the Molecule{T} associated with the given idx in sys. Throws a KeyError if no such molecule exists.

source
BiochemicalAlgorithms.moleculesFunction
molecules(::System{T} = default_system())

Returns a MoleculeTable{T} containing all molecules of the given system.

Supported keyword arguments

  • variant::Union{Nothing, MoleculeVariantType} = nothing

The keyword argument limits the results to molecules matching the given variant type. Keyword arguments set to nothing are ignored.

source
BiochemicalAlgorithms.parent_moleculeFunction
parent_molecule(::Atom)
parent_molecule(::Chain)
parent_molecule(::Fragment)

Returns the Molecule{T} containing the given object. Returns nothing if no such molecule exists.

source
BiochemicalAlgorithms.sort_molecules!Function
sort_molecules!(::System)

Sorts the molecules in the given system by idx (default) or according to the given keyword arguments.

Supported keyword arguments

Same as Base.sort

source
Base.delete!Method
delete!(::Molecule)
delete!(::MoleculeTable)
delete!(::MoleculeTable, idx::Int)

Removes the given molecule(s) and all associated chains, secondary_structures and fragments from the associated system.

Supported keyword arguments

  • keep_atoms::Bool = false Determines whether associated atoms (and their bonds) are removed as well
source
Base.push!Method
push!(::System{T}, ::Molecule{T})

Creates a copy of the given molecule in the given system. The new molecule is automatically assigned a new idx.

source

Proteins (MoleculeVariant.Protein)

Chains

BiochemicalAlgorithms.ChainType
Chain{T} <: AbstractAtomContainer{T}

Mutable representation of an individual chain in a system.

Public fields

  • idx::Int
  • name::String

Private fields

  • properties::Properties
  • flags::Flags
  • molecule_idx::Int

Constructors

Chain(
    mol::Molecule{T};
    # keyword arguments
    name::AbstractString = "",
    properties::Properties = Properties(),
    flags::Flags = Flags()
)

Creates a new Chain{T} in the given molecule.

source
BiochemicalAlgorithms.ChainTableType
ChainTable{T} <: AbstractSystemComponentTable{T}

Tables.jl-compatible representation of system chains (or a subset thereof). Chain tables can be generated using chains or filtered from other chain tables (via Base.filter).

Public columns

  • idx::AbstractVector{Int}
  • name::AbstractVector{String}

Private columns

  • properties::AbstractVector{Properties}
  • flags::AbstractVector{Flags}
  • molecule_idx::AbstractVector{Int}
source
BiochemicalAlgorithms.chain_by_idxFunction
chain_by_idx(
    sys::System{T} = default_system(),
    idx::Int
) -> Chain{T}

Returns the Chain{T} associated with the given idx in sys. Throws a KeyError if no such chain exists.

source
BiochemicalAlgorithms.chainsFunction
chains(::Molecule)
chains(::System = default_system(); kwargs...)

Returns a ChainTable{T} containing all chains of the given atom container.

Supported keyword arguments

  • molecule_idx::MaybeInt = nothing: Any value other than nothing limits the result to chains belonging to the molecule with the given ID.
source
chains(::MoleculeTable)

Returns a ChainTable{T} containing all chains of the given molecule table.

source
BiochemicalAlgorithms.nchainsFunction
nchains(::Molecule)
nchains(::System = default_system(); kwargs...)

Returns the number of chains in the given atom container.

Supported keyword arguments

See chains

source
nchains(::ChainTable)
nchains(::MolculeTable)

Returns the number of chains belonging to the given molecule table.

source
BiochemicalAlgorithms.sort_chains!Function
sort_chains!(::System)

Sorts the chains in the given system by idx (default) or according to the given keyword arguments.

Supported keyword arguments

Same as Base.sort

source
Base.delete!Method
delete!(::Chain)
delete!(::ChainTable)
delete!(::ChainTable, idx::Int)

Removes the given chain(s) and all associated secondary structures and fragments from the associated system.

Supported keyword arguments

  • keep_atoms::Bool = false Determines whether associated atoms (and their bonds) are removed as well
source
Base.push!Method
push!(::Molecule{T}, ::Chain{T})

Creates a copy of the given chain in the given molecule. The new chain is automatically assigned a new idx.

source

Secondary structures

BiochemicalAlgorithms.SecondaryStructureType
SecondaryStructure{T} <: AbstractAtomContainer{T}

Mutable representation of an individual secondary structure element in a Chain.

Public fields

  • idx::Int
  • number::Int
  • type::SecondaryStructureType
  • name::String

Private fields

  • properties::Properties
  • flags::Flags
  • molecule_idx::Int
  • chain_idx::Int

Constructors

SecondaryStructure(
    chain::Chain{T};
    number::Int,
    type::SecondaryStructureType;
    # keyword arguments
    name::AbstractString = "",
    properties::Properties = Properties(),
    flags::Flags = Flags()
)

Creates a new SecondaryStructure{T} in the given chain.

source
BiochemicalAlgorithms.SecondaryStructureTableType
SecondaryStructureTable{T} <: AbstractSystemComponentTable{T}

Tables.jl-compatible representation of system secondary structures (or a subset thereof). Secondary structure tables can be generated using secondary_structures or filtered from other secondary structure tables (via Base.filter).

Public columns

  • idx::AbstractVector{Int}
  • number::Int
  • type::SecondaryStructureType
  • name::AbstractVector{String}

Private columns

  • properties::AbstractVector{Properties}
  • flags::AbstractVector{Flags}
  • molecule_idx::AbstractVector{Int}
  • chain_idx::AbstractVector{Int}
source
BiochemicalAlgorithms.secondary_structure_by_idxFunction
secondary_structure_by_idx(
    sys::System{T},
    idx::Int64
) -> SecondaryStructure

Returns the SecondaryStructure{T} associated with the given idx in sys. Throws a KeyError if no such secondary structure exists.

source
BiochemicalAlgorithms.secondary_structuresFunction
secondary_structures(::Molecule)
secondary_structures(::Chain)
secondary_structures(::System; kwargs...)

Returns a SecondaryStructureTable{T} containing all secondary structures of the given atom container.

Supported keyword arguments

  • molecule_idx::MaybeInt = nothing
  • chain_idx::MaybeInt = nothing

All keyword arguments limit the results to secondary structures matching the given IDs. Keyword arguments set to nothing are ignored.

source
secondary_structures(::ChainTable)
secondary_structures(::MoleculeTable)

Returns a SecondaryStructureTable{T} containing all secondary structures of the given table.

Supported keyword arguments

See secondary_structures

source
BiochemicalAlgorithms.nsecondary_structuresFunction
nsecondary_structures(::Chain)
nsecondary_structures(::Molecule)
nsecondary_structures(::System; kwargs...)

Returns the number of secondary structures in the given atom container.

Supported keyword arguments

See secondary_structures

source
nsecondary_structures(::ChainTable)
nsecondary_structures(::SecondaryStructureTable)
nsecondary_structures(::MoleculeTable)

Returns the number of secondary_structures belonging to the given table.

Supported keyword arguments

See secondary_structures

source
BiochemicalAlgorithms.parent_secondary_structureFunction
parent_secondary_structure(::Atom)
parent_secondary_structure(::Fragment)
parent_secondary_structure(::Nucleotide)
parent_secondary_structure(::Residue)

Returns the SecondaryStructure{T} containing the given object. Returns nothing if no such secondary structure exists.

source
Base.delete!Method
delete!(::SecondaryStructure; keep_fragments::Bool = false)
delete!(::SecondaryStructureTable; keep_fragments::Bool = false)
delete!(::SecondaryStructureTable, idx::Int; keep_fragments::Bool = false)

Removes the given secondary_structure(s) from the associated system.

Supported keyword arguments

  • keep_fragments::Bool = false Determines whether the fragments contained in this secondary structure are removed as well. All atoms contained in those fragments are deleted as well.
source
Base.push!Method
push!(::Chain{T}, ::SecondaryStructure{T})

Creates a copy of the given secondary structure in the given chain. The new secondary structure is automatically assigned a new idx.

source

Fragments

BiochemicalAlgorithms.FragmentVariantModule
@enumx FragmentVariant begin
    None = 1
    Residue = 2
    Nucleotide = 3
end

Enum representing variants of fragments

Example

julia> res = Residue(Chain(Molecule(System())), 1; name = "ALA")
Fragment{Float32}: (idx = 3, number = 1, name = "ALA")

julia> isresidue(res)
true

julia> res.variant == FragmentVariant.Residue
true
source

Fragments (all variants)

BiochemicalAlgorithms.FragmentType
Fragment{T} <: AbstractAtomContainer{T}

Mutable representation of an individual fragment in a system.

Public fields

  • idx::Int
  • number::Int
  • name::String

Private fields

  • variant::FragmentVariantType
  • properties::Properties
  • flags::Flags
  • molecule_idx::Int
  • chain_idx::Int
  • secondary_structure_idx::Int

Constructors

Fragment(
    chain::Chain{T},
    number::Int;
    # keyword arguments
    name::AbstractString = "",
    variant::FragmentVariantType = FragmentVariant.None,
    properties::Properties = Properties(),
    flags::Flags = Flags()
)

Creates a new Fragment{T} in the given chain.

Fragment(
    secondary_structure::SecondaryStructure{T},
    number::Int;
    # keyword arguments
    name::AbstractString = "",
    variant::FragmentVariantType = FragmentVariant.None,
    properties::Properties = Properties(),
    flags::Flags = Flags()
)

Creates a new Fragment{T} in the given secondary structure.

source
BiochemicalAlgorithms.FragmentTableType
FragmentTable{T} <: AbstractSystemComponentTable{T}

Tables.jl-compatible representation of system fragments (or a subset thereof). Fragment tables can be generated using fragments or filtered from other fragment tables (via Base.filter).

Public columns

  • idx::AbstractVector{Int}
  • number::AbstractVector{Int}
  • name::AbstractVector{String}

Private columns

  • variant::AbstractVector{FragmentVariantType}
  • properties::AbstractVector{Properties}
  • flags::AbstractVector{Flags}
  • molecule_idx::AbstractVector{Int}
  • chain_idx::AbstractVector{Int}
source
BiochemicalAlgorithms.fragment_by_idxFunction
fragment_by_idx(
    sys::System{T} = default_system(),
    idx::Int
) -> Fragment{T}

Returns the Fragment{T} associated with the given idx in sys. Throws a KeyError if no such fragment exists.

source
BiochemicalAlgorithms.fragmentsFunction
fragments(::Chain)
fragments(::SecondaryStructure)
fragments(::Molecule)
fragments(::System = default_system())

Returns a FragmentTable{T} containing all fragments of the given atom container.

Supported keyword arguments

  • variant::Union{Nothing, FragmentVariantType} = nothing
  • molecule_idx::MaybeInt = nothing
  • chain_idx::MaybeInt = nothing
  • secondary_structure_idx::MaybeInt = nothing

All keyword arguments limit the results to fragments matching the given IDs or variant type. Keyword arguments set to nothing are ignored.

source
fragments(::SecondaryStructureTable)
fragments(::ChainTable)
fragments(::MoleculeTable)

Returns a FragmentTable{T} containing all fragments of the given table.

Supported keyword arguments

  • variant::Union{Nothing, FragmentVariantType} = nothing Any value other than nothing limits the results to the matching variant type.
source
BiochemicalAlgorithms.nfragmentsFunction
nfragments(::Chain)
nfragments(::SecondaryStructure)
nfragments(::Molecule)
nfragments(::System = default_system())

Returns the number of fragments in the given atom container.

Supported keyword arguments

See fragments

source
nfragments(::ChainTable)
nfragments(::FragmentTable)
nfragments(::MoleculeTable)

Returns the number of fragments belonging to the given table.

Supported keyword arguments

See fragments

source
BiochemicalAlgorithms.sort_fragments!Function
sort_fragments!(::System)

Sorts the fragments in the given system by idx (default) or according to the given keyword arguments.

Supported keyword arguments

Same as Base.sort

source
Base.delete!Method
delete!(::Fragment)
delete!(::FragmentTable)
delete!(::FragmentTable, idx::Int)

Removes the given fragment(s) from the associated system.

Supported keyword arguments

  • keep_atoms::Bool = false Determines whether associated atoms (and their bonds) are removed as well
source
Base.push!Method
push!(::SecondaryStructure{T}, ::Fragment{T})

Creates a copy of the given fragment in the given SecondaryStructure. The new fragment is automatically assigned a new idx.

source
Base.push!Method
push!(::Chain{T}, ::Fragment{T})

Creates a copy of the given fragment in the given chain. The new fragment is automatically assigned a new idx.

source

Nucleotides (FragmentVariant.Nucleotide)

BiochemicalAlgorithms.nnucleotidesFunction
nnucleotides(::Chain)
nnucleotides(::SecondaryStructure)
nnucleotides(::Molecule)
nnucleotides(::System = default_system())

Returns the number of FragmentVariant.Nucleotide fragments in the given atom container.

Supported keyword arguments

See fragments

source
nnucleotides(::ChainTable)
nnucleotides(::SecondaryStructureTable)
nnucleotides(::FragmentTable)
nnucleotides(::MoleculeTable)

Returns the number of FragmentVariant.Nucleotide fragments in the given table.

source
BiochemicalAlgorithms.nucleotidesFunction
nucleotides(::Chain)
nucleotides(::ChainTable)
nucleotides(::SecondaryStructure)
nucleotides(::SecondaryStructureTable)
nucleotides(::Molecule)
nucleotides(::MoleculeTable)
nucleotides(::System = default_system())

Returns a FragmentTable{T} containing all FragmentVariant.Nucleotide fragments of the given atom container or table.

Supported keyword arguments

See fragments

source

Nucleotides (FragmentVariant.Residue)

BiochemicalAlgorithms.residuesFunction
residues(::Chain)
residues(::ChainTable)
residues(::SecondaryStructure)
residues(::SecondaryStructureTable)
residues(::Molecule)
residues(::MoleculeTable)
residues(::System = default_system())

Returns a FragmentTable{T} containing all FragmentVariant.Residue fragments of the given atom container or table.

Supported keyword arguments

See fragments

source