Biomolecular systems
BiochemicalAlgorithms.FragmentVariantBiochemicalAlgorithms.MoleculeVariantBiochemicalAlgorithms.AbstractAtomContainerBiochemicalAlgorithms.AbstractColumnTableBiochemicalAlgorithms.AbstractSystemComponentBiochemicalAlgorithms.AbstractSystemComponentTableBiochemicalAlgorithms.AtomBiochemicalAlgorithms.AtomTableBiochemicalAlgorithms.BondBiochemicalAlgorithms.BondTableBiochemicalAlgorithms.ChainBiochemicalAlgorithms.ChainTableBiochemicalAlgorithms.FragmentBiochemicalAlgorithms.FragmentTableBiochemicalAlgorithms.FragmentVariantTypeBiochemicalAlgorithms.MoleculeBiochemicalAlgorithms.MoleculeTableBiochemicalAlgorithms.MoleculeVariantTypeBiochemicalAlgorithms.SecondaryStructureBiochemicalAlgorithms.SecondaryStructureTableBiochemicalAlgorithms.SystemBiochemicalAlgorithms.SystemComponentTableBiochemicalAlgorithms.SystemComponentTableColBase.append!Base.delete!Base.delete!Base.delete!Base.delete!Base.delete!Base.delete!Base.empty!Base.parentBase.push!Base.push!Base.push!Base.push!Base.push!Base.push!Base.push!Base.sortBase.sort!BiochemicalAlgorithms.NucleotideBiochemicalAlgorithms.ProteinBiochemicalAlgorithms.ResidueBiochemicalAlgorithms.apply_torsion_angle!BiochemicalAlgorithms.atom_by_idxBiochemicalAlgorithms.atom_by_nameBiochemicalAlgorithms.atomsBiochemicalAlgorithms.bond_by_idxBiochemicalAlgorithms.bondsBiochemicalAlgorithms.calculate_bond_angleBiochemicalAlgorithms.calculate_torsion_angleBiochemicalAlgorithms.chain_by_idxBiochemicalAlgorithms.chainsBiochemicalAlgorithms.default_systemBiochemicalAlgorithms.fragment_by_idxBiochemicalAlgorithms.fragmentsBiochemicalAlgorithms.full_tableBiochemicalAlgorithms.get_propertyBiochemicalAlgorithms.get_torsion_omegaBiochemicalAlgorithms.get_torsion_phiBiochemicalAlgorithms.get_torsion_psiBiochemicalAlgorithms.has_flagBiochemicalAlgorithms.has_propertyBiochemicalAlgorithms.has_torsion_omegaBiochemicalAlgorithms.has_torsion_phiBiochemicalAlgorithms.has_torsion_psiBiochemicalAlgorithms.is_bound_toBiochemicalAlgorithms.is_geminalBiochemicalAlgorithms.is_vicinalBiochemicalAlgorithms.isnucleotideBiochemicalAlgorithms.isproteinBiochemicalAlgorithms.isresidueBiochemicalAlgorithms.molecule_by_idxBiochemicalAlgorithms.moleculesBiochemicalAlgorithms.natomsBiochemicalAlgorithms.nbondsBiochemicalAlgorithms.nchainsBiochemicalAlgorithms.nfragmentsBiochemicalAlgorithms.nmoleculesBiochemicalAlgorithms.nnucleotidesBiochemicalAlgorithms.nproteinsBiochemicalAlgorithms.nresiduesBiochemicalAlgorithms.nsecondary_structuresBiochemicalAlgorithms.nucleotide_by_idxBiochemicalAlgorithms.nucleotidesBiochemicalAlgorithms.parent_chainBiochemicalAlgorithms.parent_fragmentBiochemicalAlgorithms.parent_moleculeBiochemicalAlgorithms.parent_nucleotideBiochemicalAlgorithms.parent_proteinBiochemicalAlgorithms.parent_residueBiochemicalAlgorithms.parent_secondary_structureBiochemicalAlgorithms.parent_systemBiochemicalAlgorithms.protein_by_idxBiochemicalAlgorithms.proteinsBiochemicalAlgorithms.residue_by_idxBiochemicalAlgorithms.residuesBiochemicalAlgorithms.revalidate_indices!BiochemicalAlgorithms.secondary_structure_by_idxBiochemicalAlgorithms.secondary_structuresBiochemicalAlgorithms.set_flag!BiochemicalAlgorithms.set_property!BiochemicalAlgorithms.sort_atoms!BiochemicalAlgorithms.sort_bonds!BiochemicalAlgorithms.sort_chains!BiochemicalAlgorithms.sort_fragments!BiochemicalAlgorithms.sort_molecules!BiochemicalAlgorithms.sort_secondary_structures!BiochemicalAlgorithms.unset_flag!
Systems
BiochemicalAlgorithms.System — Type
mutable struct System{T} <: AbstractAtomContainer{T}Mutable representation of a biomolecular system.
Fields
name::Stringproperties::Propertiesflags::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}.
BiochemicalAlgorithms.default_system — Function
default_system() -> System{Float32}
Returns the global default system.
BiochemicalAlgorithms.parent_system — Function
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.
Base.parent — Method
parent(::Atom)
parent(::Bond)
parent(::Chain)
parent(::SecondaryStructure)
parent(::Fragment)
parent(::Molecule)
parent(::System)Returns the System{T} containing the given object.
Base.append! — Method
Base.append!(sys::System, others::System...)Copies all system components from others to sys.
Base.empty! — Method
empty!(::System)Removes all components from the system.
System components
BiochemicalAlgorithms.AbstractAtomContainer — Type
abstract type AbstractAtomContainer{T} <: AbstractSystemComponent{T}Abstract base type for all atom containers.
BiochemicalAlgorithms.AbstractColumnTable — Type
abstract type AbstractColumnTable <: Tables.AbstractColumnsAbstract base type for all Tables.jl-compatible column tables.
BiochemicalAlgorithms.AbstractSystemComponent — Type
abstract type AbstractSystemComponent{T<:Real}Abstract base type for all components of a system, including the system itself.
BiochemicalAlgorithms.AbstractSystemComponentTable — Type
abstract type AbstractSystemComponentTable{T<:Real} <: AbstractColumnTableAbstract base type for all Tables.jl-compatible system component tables.
BiochemicalAlgorithms.SystemComponentTable — Type
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.).
BiochemicalAlgorithms.SystemComponentTableCol — Type
struct SystemComponentTableCol{T} <: AbstractArray{T, 1}Vector-like representation of a single SystemComponentTable column.
BiochemicalAlgorithms.full_table — Function
full_table(::SystemComponentTable)Returns an extended copy of the given table, with all columns being visible.
BiochemicalAlgorithms.get_property — Function
get_property(
ac::AbstractSystemComponent,
key::Symbol
) -> Any
Returns the property associated with the given key in ac.
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.
BiochemicalAlgorithms.has_flag — Function
has_flag(ac::AbstractSystemComponent, flag::Symbol) -> Any
Returns a Bool indicating whether the given system component has the given flag.
BiochemicalAlgorithms.has_property — Function
has_property(
ac::AbstractSystemComponent,
key::Symbol
) -> Any
Returns a Bool indicating whether the given system component has the given property.
BiochemicalAlgorithms.revalidate_indices! — Function
revalidate_indices!(::SystemComponentTable)
revalidate_indices!(::SystemComponentTableCol)Removes remnants of previously removed system components from the given table or table column.
BiochemicalAlgorithms.set_flag! — Function
set_flag!(ac::AbstractSystemComponent, flag::Symbol) -> Any
Adds the given flag to ac.
BiochemicalAlgorithms.set_property! — Function
set_property!(
ac::AbstractSystemComponent,
key::Symbol,
value
) -> Any
Sets the property associated with the given key in ac to the given value.
BiochemicalAlgorithms.unset_flag! — Function
unset_flag!(
ac::AbstractSystemComponent,
flag::Symbol
) -> Any
Removes the given flag from ac.
Base.sort! — Method
sort!(::SystemComponentTable)In-place variant of sort.
Atoms
BiochemicalAlgorithms.Atom — Type
Atom{T} <: AbstractSystemComponent{T}Mutable representation of an individual atom in a system.
Public fields
idx::Intnumber::Intelement::ElementTypename::Stringatom_type::Stringr::Vector3{T}v::Vector3{T}F::Vector3{T}formal_charge::Intcharge::Tradius::T
Private fields
properties::Propertiesflags::Flagsframe_id::Intmolecule_idx::MaybeIntchain_idx::MaybeIntfragment_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.
BiochemicalAlgorithms.AtomTable — Type
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}
BiochemicalAlgorithms.atom_by_idx — Function
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.
BiochemicalAlgorithms.atom_by_name — Function
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 thannothinglimits the result to atoms matching this frame ID.
BiochemicalAlgorithms.atoms — Function
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 = 1molecule_idx::Union{MaybeInt, Some{Nothing}} = nothingchain_idx::Union{MaybeInt, Some{Nothing}} = nothingfragment_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.
atoms(::ChainTable)
atoms(::FragmentTable)
atoms(::MoleculeTable)Returns an AtomTable{T} containing all atoms of the given table.
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.
atoms(
ab::AbstractAtomBijection
) -> Tuple{AtomTable, AtomTable}
Returns the tuple of atom tables represented by this bijection.
BiochemicalAlgorithms.is_bound_to — Function
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.
BiochemicalAlgorithms.is_geminal — Function
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.
BiochemicalAlgorithms.is_vicinal — Function
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.
BiochemicalAlgorithms.natoms — Function
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
natoms(::AtomTable)
natoms(::ChainTable)
natoms(::FragmentTable)
natoms(::MoleculeTable)Returns the number of atoms in the given table.
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
Base.delete! — Method
delete!(::Atom)
delete!(::AtomTable)
delete!(::AtomTable, idx::Int)Removes the given atom(s) and all associated bonds from the associated system.
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
Bonds
BiochemicalAlgorithms.Bond — Type
Bond{T} <: AbstractAtomContainer{T}Mutable representation of an individual bond in a system.
Public fields
idx::Inta1::Inta2::Intorder::BondOrderType
Private fields
properties::Propertiesflags::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.
BiochemicalAlgorithms.BondTable — Type
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}
BiochemicalAlgorithms.bond_by_idx — Function
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.
BiochemicalAlgorithms.bonds — Function
bonds(::Atom)
bonds(::AtomTable)Returns a BondTable{T} containing all bonds of the given atom(s).
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
bonds(::ChainTable)
bonds(::FragmentTable)
bonds(::MoleculeTable)Returns a BondTable{T} containing all bonds where at least one associated atom belongs to the given table.
BiochemicalAlgorithms.nbonds — Function
nbonds(::Atom)
nbonds(::AtomTable)Returns the number of bonds of the given atom(s).
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
nbonds(::BondTable)
nbonds(::ChainTable)
nbonds(::FragmentTable)
nbonds(::MoleculeTable)Returns the number of bonds where at least one associated atom belongs to the given table.
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
Base.delete! — Method
delete!(::Bond)
delete!(::BondTable)
delete!(::BondTable, idx::Int)Removes the given bond(s) from the associated system.
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.
Molecules
BiochemicalAlgorithms.MoleculeVariant — Module
@enumx MoleculeVariant begin
None = 1
Protein = 2
endEnum representing variants of molecules
Example
julia> prot = Protein(System())
Molecule{Float32}: (idx = 1, name = "")
julia> isprotein(prot)
true
julia> prot.variant == MoleculeVariant.Protein
trueBiochemicalAlgorithms.MoleculeVariantType — Type
const MoleculeVariantType = MoleculeVariant.TType of MoleculeVariant enumerators
Molecules (all variants)
BiochemicalAlgorithms.Molecule — Type
Molecule{T} <: AbstractAtomContainer{T}Mutable representation of an individual molecule in a system.
Public fields
idx::Intname::String
Private fields
variant::MoleculeVariantTypeproperties::Propertiesflags::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.
BiochemicalAlgorithms.MoleculeTable — Type
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}
BiochemicalAlgorithms.molecule_by_idx — Function
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.
BiochemicalAlgorithms.molecules — Function
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.
BiochemicalAlgorithms.nmolecules — Function
nmolecules(::System = default_system())Returns the number of molecules in the given system.
Supported keyword arguments
See molecules
nmolecules(::MoleculeTable)Returns the number of molecules in the given table.
Supported keyword arguments
See molecules
BiochemicalAlgorithms.parent_molecule — Function
parent_molecule(::Atom)
parent_molecule(::Chain)
parent_molecule(::Fragment)Returns the Molecule{T} containing the given object. Returns nothing if no such molecule exists.
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
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 = falseDetermines whether associated atoms (and their bonds) are removed as well
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.
Proteins (MoleculeVariant.Protein)
BiochemicalAlgorithms.Protein — Function
Protein(sys::System = default_system())Molecule constructor defaulting to the MoleculeVariant.Protein variant.
Supported keyword arguments
See Molecule
BiochemicalAlgorithms.isprotein — Function
isprotein(mol::Molecule) -> Bool
Returns true if the given molecule is a MoleculeVariant.Protein, false otherwise.
BiochemicalAlgorithms.nproteins — Function
nproteins(::System = default_system())Returns the number of MoleculeVariant.Protein molecules in the given system.
Supported keyword arguments
See molecules
nproteins(::MoleculeTable)Returns the number of MoleculeVariant.Protein molecules in the given table.
BiochemicalAlgorithms.parent_protein — Function
parent_protein(::Atom)
parent_protein(::Chain)
parent_protein(::Fragment)Returns the Molecule{T} containing the given object. Returns nothing if no such molecule exists or if the molecule is not a MoleculeVariant.Protein.
BiochemicalAlgorithms.protein_by_idx — Function
protein_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 or if the molecule is not a MoleculeVariant.Protein.
BiochemicalAlgorithms.proteins — Function
proteins(::System{T} = default_system())Returns a MoleculeTable{T} containing all MoleculeVariant.Protein molecules of the given system.
Supported keyword arguments
See molecules
Chains
BiochemicalAlgorithms.Chain — Type
Chain{T} <: AbstractAtomContainer{T}Mutable representation of an individual chain in a system.
Public fields
idx::Intname::String
Private fields
properties::Propertiesflags::Flagsmolecule_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.
BiochemicalAlgorithms.ChainTable — Type
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}
BiochemicalAlgorithms.chain_by_idx — Function
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.
BiochemicalAlgorithms.chains — Function
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 thannothinglimits the result to chains belonging to the molecule with the given ID.
chains(::MoleculeTable)Returns a ChainTable{T} containing all chains of the given molecule table.
BiochemicalAlgorithms.nchains — Function
nchains(::Molecule)
nchains(::System = default_system(); kwargs...)Returns the number of chains in the given atom container.
Supported keyword arguments
See chains
nchains(::ChainTable)
nchains(::MolculeTable)Returns the number of chains belonging to the given molecule table.
BiochemicalAlgorithms.parent_chain — Function
parent_chain(::Atom)
parent_chain(::Fragment)Returns the Chain{T} containing the given object. Returns nothing if no such chain exists.
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
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 = falseDetermines whether associated atoms (and their bonds) are removed as well
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.
Secondary structures
BiochemicalAlgorithms.SecondaryStructure — Type
SecondaryStructure{T} <: AbstractAtomContainer{T}Mutable representation of an individual secondary structure element in a Chain.
Public fields
idx::Intnumber::Inttype::SecondaryStructureTypename::String
Private fields
properties::Propertiesflags::Flagsmolecule_idx::Intchain_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.
BiochemicalAlgorithms.SecondaryStructureTable — Type
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::Inttype::SecondaryStructureTypename::AbstractVector{String}
Private columns
properties::AbstractVector{Properties}flags::AbstractVector{Flags}molecule_idx::AbstractVector{Int}chain_idx::AbstractVector{Int}
BiochemicalAlgorithms.secondary_structure_by_idx — Function
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.
BiochemicalAlgorithms.secondary_structures — Function
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 = nothingchain_idx::MaybeInt = nothing
All keyword arguments limit the results to secondary structures matching the given IDs. Keyword arguments set to nothing are ignored.
secondary_structures(::ChainTable)
secondary_structures(::MoleculeTable)Returns a SecondaryStructureTable{T} containing all secondary structures of the given table.
Supported keyword arguments
BiochemicalAlgorithms.nsecondary_structures — Function
nsecondary_structures(::Chain)
nsecondary_structures(::Molecule)
nsecondary_structures(::System; kwargs...)Returns the number of secondary structures in the given atom container.
Supported keyword arguments
nsecondary_structures(::ChainTable)
nsecondary_structures(::SecondaryStructureTable)
nsecondary_structures(::MoleculeTable)Returns the number of secondary_structures belonging to the given table.
Supported keyword arguments
BiochemicalAlgorithms.parent_secondary_structure — Function
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.
BiochemicalAlgorithms.sort_secondary_structures! — Function
sort_secondary_structures!(::System)Sorts the secondary structures in the given system by idx (default) or according to the given keyword arguments.
Supported keyword arguments
Same as Base.sort
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 = falseDetermines whether the fragments contained in this secondary structure are removed as well. All atoms contained in those fragments are deleted as well.
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.
Fragments
BiochemicalAlgorithms.FragmentVariant — Module
@enumx FragmentVariant begin
None = 1
Residue = 2
Nucleotide = 3
endEnum 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
trueBiochemicalAlgorithms.FragmentVariantType — Type
const FragmentVariantType = FragmentVariant.TType of FragmentVariant enumerators
BiochemicalAlgorithms.apply_torsion_angle! — Function
apply_torsion_angle!(a::Atom{T}, b::Atom{T}, c::Atom{T}, d::Atom{T}, angle::T)Applies the torsion angle defined by the four atoms a, b, c, and d to the specified angle in radians. Returns true if the torsion angle has been applied successfully, false otherwise.
This function modifies the positions of the atoms to achieve the specified torsion angle. It assumes that the atoms are part of a molecular structure where such modifications are meaningful. It does not check for steric clashes or other geometric or chemical constraints that may arise from the rotation.
BiochemicalAlgorithms.calculate_bond_angle — Function
calculate_bond_angle(a::Atom{T}, b::Atom{T}, c::Atom{T}) -> TCalculates the bond angle formed by three atoms a, b, and c and returns it (in radians). The angle is measured at atom b, with a and c being the other two atoms forming the angle.
BiochemicalAlgorithms.calculate_torsion_angle — Function
calculate_torsion_angle(a::Atom{T}, b::Atom{T}, c::Atom{T}, d::Atom{T}) -> TCalculates the torsion angle (dihedral angle, in radians) defined by four atoms a, b, c, and d. The torsion angle is the angle between the plane formed by atoms a, b, c and the plane formed by atoms b, c, d.
BiochemicalAlgorithms.get_torsion_omega — Function
get_torsion_omega(frag::Fragment{T})) -> TCalculates the torsion angle omega (in radians) defined by the given frag::Fragment object. Returns the calculated torsion angle, if the angle can be calculated, else zero. The omega angle is defined by the fragment's CA and C atoms, and the next fragment's CA and N atoms.
BiochemicalAlgorithms.get_torsion_phi — Function
get_torsion_phi(frag::Fragment{T}) -> TCalculates the torsion angle phi (in radians) defined by the given frag::Fragment object. Returns the calculated torsion angle, if the angle can be calculated, else zero. The phi angle is defined by the fragment's atoms N, CA, C, and the previous fragment's C atom.
BiochemicalAlgorithms.get_torsion_psi — Function
get_torsion_psi(frag::Fragment{T}) -> TCalculates the torsion angle psi (in radians) defined by the given frag::Fragment object. Returns the calculated torsion angle, if the angle can be calculated, else zero. The psi angle is defined by the fragment's atoms N, CA, C, and the next fragment's N atom.
BiochemicalAlgorithms.has_torsion_omega — Function
has_torsion_omega(frag::Fragment) -> BoolChecks if the given Fragment object has a torsion angle omega defined, false otherwise. See Wikipedia for more information.
BiochemicalAlgorithms.has_torsion_phi — Function
has_torsion_phi(frag::Fragment) -> BoolChecks if the given frag::Fragment object has a torsion angle phi defined, false otherwise. See Wikipedia for more information.
BiochemicalAlgorithms.has_torsion_psi — Function
has_torsion_psi(frag::Fragment) -> BoolChecks if the given frag::Fragment object has a torsion angle psi defined, false otherwise. See Wikipedia for more information.
Fragments (all variants)
BiochemicalAlgorithms.Fragment — Type
Fragment{T} <: AbstractAtomContainer{T}Mutable representation of an individual fragment in a system.
Public fields
idx::Intnumber::Intname::String
Private fields
variant::FragmentVariantTypeproperties::Propertiesflags::Flagsmolecule_idx::Intchain_idx::Intsecondary_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.
BiochemicalAlgorithms.FragmentTable — Type
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}
BiochemicalAlgorithms.fragment_by_idx — Function
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.
BiochemicalAlgorithms.fragments — Function
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} = nothingmolecule_idx::MaybeInt = nothingchain_idx::MaybeInt = nothingsecondary_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.
fragments(::SecondaryStructureTable)
fragments(::ChainTable)
fragments(::MoleculeTable)Returns a FragmentTable{T} containing all fragments of the given table.
Supported keyword arguments
variant::Union{Nothing, FragmentVariantType} = nothingAny value other thannothinglimits the results to the matching variant type.
BiochemicalAlgorithms.nfragments — Function
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
nfragments(::ChainTable)
nfragments(::FragmentTable)
nfragments(::MoleculeTable)Returns the number of fragments belonging to the given table.
Supported keyword arguments
See fragments
BiochemicalAlgorithms.parent_fragment — Function
parent_fragment(::Atom)Returns the Fragment{T} containing the given atom. Returns nothing if no such fragment exists.
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
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 = falseDetermines whether associated atoms (and their bonds) are removed as well
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.
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.
Nucleotides (FragmentVariant.Nucleotide)
BiochemicalAlgorithms.Nucleotide — Function
Nucleotide(chain::Chain, number::Int)Fragment constructor defaulting to the FragmentVariant.Nucleotide variant.
Supported keyword arguments
See Fragment
Nucleotide(ss::SecondaryStructure, number::Int)Fragment constructor defaulting to the FragmentVariant.Nucleotide variant.
Supported keyword arguments
See Fragment
BiochemicalAlgorithms.isnucleotide — Function
isnucleotide(frag::Fragment) -> Bool
Returns true if the given fragment is a FragmentVariant.Nucleotide, false otherwise.
BiochemicalAlgorithms.nnucleotides — Function
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
nnucleotides(::ChainTable)
nnucleotides(::SecondaryStructureTable)
nnucleotides(::FragmentTable)
nnucleotides(::MoleculeTable)Returns the number of FragmentVariant.Nucleotide fragments in the given table.
BiochemicalAlgorithms.nucleotide_by_idx — Function
nucleotide_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 or if the fragment is not a FragmentVariant.Nucleotide.
BiochemicalAlgorithms.nucleotides — Function
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
BiochemicalAlgorithms.parent_nucleotide — Function
parent_nucleotide(::Atom)Returns the Fragment{T} containing the given atom. Returns nothing if no such fragment exists or if the fragment is not a FragmentVariant.Nucleotide.
Nucleotides (FragmentVariant.Residue)
BiochemicalAlgorithms.Residue — Function
Residue(chain::Chain, number::Int)Fragment constructor defaulting to the FragmentVariant.Residue variant.
Supported keyword arguments
See Fragment
Residue(ss::SecondaryStructure, number::Int)Fragment constructor defaulting to the FragmentVariant.Residue variant.
Supported keyword arguments
See Fragment
BiochemicalAlgorithms.isresidue — Function
isresidue(frag::Fragment) -> Bool
Returns true if the given fragment is a FragmentVariant.Residue, false otherwise.
BiochemicalAlgorithms.nresidues — Function
nresidues(::Chain)
nresidues(::SecondaryStructure)
nresidues(::Molecule)
nresidues(::System = default_system())Returns the number of FragmentVariant.Residue fragments in the given atom container.
Supported keyword arguments
See fragments
nresidues(::ChainTable)
nresidues(::SecondaryStructureTable)
nresidues(::FragmentTable)
nresidues(::MoleculeTable)Returns the number of FragmentVariant.Residue fragments in the given table.
BiochemicalAlgorithms.parent_residue — Function
parent_residue(::Atom)Returns the Fragment{T} containing the given atom. Returns nothing if no such fragment exists or if the fragment is not a FragmentVariant.Residue.
BiochemicalAlgorithms.residue_by_idx — Function
residue_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 or if the fragment is not a FragmentVariant.Residue.
BiochemicalAlgorithms.residues — Function
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