Biomolecular systems
BiochemicalAlgorithms.FragmentVariant
BiochemicalAlgorithms.MoleculeVariant
BiochemicalAlgorithms.AbstractAtomContainer
BiochemicalAlgorithms.AbstractColumnTable
BiochemicalAlgorithms.AbstractSystemComponent
BiochemicalAlgorithms.AbstractSystemComponentTable
BiochemicalAlgorithms.Atom
BiochemicalAlgorithms.AtomTable
BiochemicalAlgorithms.Bond
BiochemicalAlgorithms.BondTable
BiochemicalAlgorithms.Chain
BiochemicalAlgorithms.ChainTable
BiochemicalAlgorithms.Fragment
BiochemicalAlgorithms.FragmentTable
BiochemicalAlgorithms.FragmentVariantType
BiochemicalAlgorithms.Molecule
BiochemicalAlgorithms.MoleculeTable
BiochemicalAlgorithms.MoleculeVariantType
BiochemicalAlgorithms.SecondaryStructure
BiochemicalAlgorithms.SecondaryStructureTable
BiochemicalAlgorithms.System
BiochemicalAlgorithms.SystemComponentTable
BiochemicalAlgorithms.SystemComponentTableCol
Base.append!
Base.delete!
Base.delete!
Base.delete!
Base.delete!
Base.delete!
Base.delete!
Base.empty!
Base.parent
Base.push!
Base.push!
Base.push!
Base.push!
Base.push!
Base.push!
Base.push!
Base.sort
Base.sort!
BiochemicalAlgorithms.Nucleotide
BiochemicalAlgorithms.Protein
BiochemicalAlgorithms.Residue
BiochemicalAlgorithms.atom_by_idx
BiochemicalAlgorithms.atom_by_name
BiochemicalAlgorithms.atoms
BiochemicalAlgorithms.bond_by_idx
BiochemicalAlgorithms.bonds
BiochemicalAlgorithms.chain_by_idx
BiochemicalAlgorithms.chains
BiochemicalAlgorithms.default_system
BiochemicalAlgorithms.fragment_by_idx
BiochemicalAlgorithms.fragments
BiochemicalAlgorithms.full_table
BiochemicalAlgorithms.get_property
BiochemicalAlgorithms.has_flag
BiochemicalAlgorithms.has_property
BiochemicalAlgorithms.is_bound_to
BiochemicalAlgorithms.is_geminal
BiochemicalAlgorithms.is_vicinal
BiochemicalAlgorithms.isnucleotide
BiochemicalAlgorithms.isprotein
BiochemicalAlgorithms.isresidue
BiochemicalAlgorithms.molecule_by_idx
BiochemicalAlgorithms.molecules
BiochemicalAlgorithms.natoms
BiochemicalAlgorithms.nbonds
BiochemicalAlgorithms.nchains
BiochemicalAlgorithms.nfragments
BiochemicalAlgorithms.nmolecules
BiochemicalAlgorithms.nnucleotides
BiochemicalAlgorithms.nproteins
BiochemicalAlgorithms.nresidues
BiochemicalAlgorithms.nsecondary_structures
BiochemicalAlgorithms.nucleotide_by_idx
BiochemicalAlgorithms.nucleotides
BiochemicalAlgorithms.parent_chain
BiochemicalAlgorithms.parent_fragment
BiochemicalAlgorithms.parent_molecule
BiochemicalAlgorithms.parent_nucleotide
BiochemicalAlgorithms.parent_protein
BiochemicalAlgorithms.parent_residue
BiochemicalAlgorithms.parent_secondary_structure
BiochemicalAlgorithms.parent_system
BiochemicalAlgorithms.protein_by_idx
BiochemicalAlgorithms.proteins
BiochemicalAlgorithms.residue_by_idx
BiochemicalAlgorithms.residues
BiochemicalAlgorithms.revalidate_indices!
BiochemicalAlgorithms.secondary_structure_by_idx
BiochemicalAlgorithms.secondary_structures
BiochemicalAlgorithms.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
— Typemutable 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}
.
BiochemicalAlgorithms.default_system
— Functiondefault_system() -> System{Float32}
Returns the global default system.
BiochemicalAlgorithms.parent_system
— Functionparent_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
— Methodparent(::Atom)
parent(::Bond)
parent(::Chain)
parent(::SecondaryStructure)
parent(::Fragment)
parent(::Molecule)
parent(::System)
Returns the System{T}
containing the given object.
Base.append!
— MethodBase.append!(sys::System, others::System...)
Copies all system components from others
to sys
.
All copied components are assigned new idx
values. The (system) names, (system) properties and (system) flags of others
are ignored.
Base.empty!
— Methodempty!(::System)
Removes all components from the system.
System components
BiochemicalAlgorithms.AbstractAtomContainer
— Typeabstract type AbstractAtomContainer{T} <: AbstractSystemComponent{T}
Abstract base type for all atom containers.
BiochemicalAlgorithms.AbstractColumnTable
— Typeabstract type AbstractColumnTable <: Tables.AbstractColumns
Abstract base type for all Tables.jl-compatible column tables.
BiochemicalAlgorithms.AbstractSystemComponent
— Typeabstract type AbstractSystemComponent{T<:Real}
Abstract base type for all components of a system, including the system itself.
BiochemicalAlgorithms.AbstractSystemComponentTable
— Typeabstract type AbstractSystemComponentTable{T<:Real} <: AbstractColumnTable
Abstract base type for all Tables.jl-compatible system component tables.
BiochemicalAlgorithms.SystemComponentTable
— Typestruct 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
— Typestruct SystemComponentTableCol{T} <: AbstractArray{T, 1}
Vector
-like representation of a single SystemComponentTable
column.
BiochemicalAlgorithms.full_table
— Functionfull_table(::SystemComponentTable)
Returns an extended copy of the given table, with all columns being visible.
BiochemicalAlgorithms.get_property
— Functionget_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
— Functionhas_flag(ac::AbstractSystemComponent, flag::Symbol) -> Any
Returns a Bool
indicating whether the given system component has the given flag.
BiochemicalAlgorithms.has_property
— Functionhas_property(
ac::AbstractSystemComponent,
key::Symbol
) -> Any
Returns a Bool
indicating whether the given system component has the given property.
BiochemicalAlgorithms.revalidate_indices!
— Functionrevalidate_indices!(::SystemComponentTable)
revalidate_indices!(::SystemComponentTableCol)
Removes remnants of previously removed system components from the given table or table column.
BiochemicalAlgorithms.set_flag!
— Functionset_flag!(ac::AbstractSystemComponent, flag::Symbol) -> Any
Adds the given flag to ac
.
BiochemicalAlgorithms.set_property!
— Functionset_property!(
ac::AbstractSystemComponent,
key::Symbol,
value
) -> Any
Sets the property associated with the given key in ac
to the given value
.
BiochemicalAlgorithms.unset_flag!
— Functionunset_flag!(
ac::AbstractSystemComponent,
flag::Symbol
) -> Any
Removes the given flag from ac
.
Base.sort
— Methodsort(::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
Base.sort!
— Methodsort!(::SystemComponentTable)
In-place variant of sort
.
Only the given table is modified, not the underlying system!
Atoms
BiochemicalAlgorithms.Atom
— TypeAtom{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.
BiochemicalAlgorithms.AtomTable
— TypeAtomTable{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
— Functionatom_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
— Functionatom_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 thannothing
limits the result to atoms matching this frame ID.
BiochemicalAlgorithms.atoms
— Functionatoms(::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
.
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
— Functionis_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
— Functionis_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
— Functionis_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
— Functionnatoms(::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!
— Functionsort_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!
— Methoddelete!(::Atom)
delete!(::AtomTable)
delete!(::AtomTable, idx::Int)
Removes the given atom(s) and all associated bonds from the associated system.
Base.push!
— Methodpush!(::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
— TypeBond{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.
BiochemicalAlgorithms.BondTable
— TypeBondTable{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
— Functionbond_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
— Functionbonds(::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
— Functionnbonds(::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!
— Functionsort_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!
— Methoddelete!(::Bond)
delete!(::BondTable)
delete!(::BondTable, idx::Int)
Removes the given bond(s) from the associated system.
Base.push!
— Methodpush!(::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
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
BiochemicalAlgorithms.MoleculeVariantType
— Typeconst MoleculeVariantType = MoleculeVariant.T
Type of MoleculeVariant
enumerators
Molecules (all variants)
BiochemicalAlgorithms.Molecule
— TypeMolecule{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.
BiochemicalAlgorithms.MoleculeTable
— TypeMoleculeTable{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
— Functionmolecule_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
— Functionmolecules(::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
— Functionnmolecules(::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
— Functionparent_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!
— Functionsort_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!
— Methoddelete!(::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
Base.push!
— Methodpush!(::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
— FunctionProtein(sys::System = default_system())
Molecule
constructor defaulting to the MoleculeVariant.Protein
variant.
Supported keyword arguments
See Molecule
BiochemicalAlgorithms.isprotein
— Functionisprotein(mol::Molecule) -> Bool
Returns true
if the given molecule is a MoleculeVariant.Protein
, false
otherwise.
BiochemicalAlgorithms.nproteins
— Functionnproteins(::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
— Functionparent_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
— Functionprotein_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
— Functionproteins(::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
— TypeChain{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.
BiochemicalAlgorithms.ChainTable
— TypeChainTable{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
— Functionchain_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
— Functionchains(::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 thannothing
limits 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
— Functionnchains(::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
— Functionparent_chain(::Atom)
parent_chain(::Fragment)
Returns the Chain{T}
containing the given object. Returns nothing
if no such chain exists.
BiochemicalAlgorithms.sort_chains!
— Functionsort_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!
— Methoddelete!(::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
Base.push!
— Methodpush!(::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
— TypeSecondaryStructure{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.
BiochemicalAlgorithms.SecondaryStructureTable
— TypeSecondaryStructureTable{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}
BiochemicalAlgorithms.secondary_structure_by_idx
— Functionsecondary_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
— Functionsecondary_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.
secondary_structures(::ChainTable)
secondary_structures(::MoleculeTable)
Returns a SecondaryStructureTable{T}
containing all secondary structures of the given table.
Supported keyword arguments
BiochemicalAlgorithms.nsecondary_structures
— Functionnsecondary_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
— Functionparent_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!
— Functionsort_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!
— Methoddelete!(::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.
Base.push!
— Methodpush!(::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
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
BiochemicalAlgorithms.FragmentVariantType
— Typeconst FragmentVariantType = FragmentVariant.T
Type of FragmentVariant
enumerators
Fragments (all variants)
BiochemicalAlgorithms.Fragment
— TypeFragment{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.
BiochemicalAlgorithms.FragmentTable
— TypeFragmentTable{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
— Functionfragment_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
— Functionfragments(::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.
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 thannothing
limits the results to the matching variant type.
BiochemicalAlgorithms.nfragments
— Functionnfragments(::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
— Functionparent_fragment(::Atom)
Returns the Fragment{T}
containing the given atom. Returns nothing
if no such fragment exists.
BiochemicalAlgorithms.sort_fragments!
— Functionsort_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!
— Methoddelete!(::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
Base.push!
— Methodpush!(::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!
— Methodpush!(::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
— FunctionNucleotide(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
— Functionisnucleotide(frag::Fragment) -> Bool
Returns true
if the given fragment is a FragmentVariant.Nucleotide
, false
otherwise.
BiochemicalAlgorithms.nnucleotides
— Functionnnucleotides(::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
— Functionnucleotide_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
— Functionnucleotides(::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
— Functionparent_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
— FunctionResidue(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
— Functionisresidue(frag::Fragment) -> Bool
Returns true
if the given fragment is a FragmentVariant.Residue
, false
otherwise.
BiochemicalAlgorithms.nresidues
— Functionnresidues(::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
— Functionparent_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
— Functionresidue_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
— Functionresidues(::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