MolSimToolkitShared.jl
MolSimToolkitShared.jl is a low-level package defining a series of functions and function names that are shared among other packages, such as MolSimToolkit.jl, PDBTools.jl, and ComplexMixtures.jl.
Normally this package won't be used by the end-user, but as a dependency of other packages.
API
Function name placeholders
These function names are considered public:
distance, distances
coordination_number, bulk_coordination
center_of_mass
wrap, wrap_to_first
align, align!, rmsd
Methods
MolSimToolkitShared.align
— Functionalign(x, y; mass = nothing)
align!(x, y; mass = nothing)
Aligns two structures (sets of points in 3D space). Solves the "Procrustes" problem, which is to find the best translation, and rotation, that aligns the two structures, minimizing the RMSD between them.
Structures are expected to be of the same size, and the correspondence is assumed from the vector indices.
align
returns a new vector containing the coordinates of x aligned to y. align!
modifies the input vector x
in place.
MolSimToolkitShared.align!
— Functionalign(x, y; mass = nothing)
align!(x, y; mass = nothing)
Aligns two structures (sets of points in 3D space). Solves the "Procrustes" problem, which is to find the best translation, and rotation, that aligns the two structures, minimizing the RMSD between them.
Structures are expected to be of the same size, and the correspondence is assumed from the vector indices.
align
returns a new vector containing the coordinates of x aligned to y. align!
modifies the input vector x
in place.
MolSimToolkitShared.center_of_mass
— Methodcenter_of_mass(x::AbstractVector{<:AbstractVector}[, mass::AbstractVector=nothing])
Calculate the center of mass of a set of points.
Arguments
x::AbstractVector{<:AbstractVector}
: A vector of coordinates.mass::AbstractVector
: A vector of masses. If not provided, all masses are assumed to be equal.
Example
julia> import MolSimToolkitShared: center_of_mass
julia> x = [ [1.0, 2.0, 3.0], [4.0, 5.0, 6.0], [7.0, 8.0, 9.0] ];
julia> center_of_mass(x)
3-element Vector{Float64}:
4.0
5.0
6.0
julia> center_of_mass(x, [1.0, 2.0, 3.0]) # providing masses
3-element Vector{Float64}:
5.0
6.0
7.0
MolSimToolkitShared.rmsd
— Methodrmsd(x::AbstractVector,y::AbstractVector)
Calculate the root mean square deviation between two vectors of coordinates.
Arguments
x::AbstractVector
: A vector of coordinates.y::AbstractVector
: A vector of coordinates.
Returns
rmsd::Real
: The root mean square deviation between the two vectors.
julia> import MolSimToolkitShared: rmsd
julia> x = [ [1.0, 2.0, 3.0], [4.0, 5.0, 6.0], [7.0, 8.0, 9.0] ];
julia> y = [ [2.0, 3.0, 4.0], [4.0, 5.0, 6.0], [7.0, 8.0, 9.0] ];
julia> rmsd(x, y)
1.0
MolSimToolkitShared.wrap
— Functionwrap(x, xref, unit_cell_matrix::SMatrix{N,N,T}) where {N,T}
wrap(x, xref, sides::AbstractVector)
Wraps the coordinates of point x
such that it is the minimum image relative to xref
. The unit cell may be given a a static matrix of size (N,N)
or as a vector of length N
.
MolSimToolkitShared.wrap_to_first
— Functionwrap_to_first(x, unit_cell_matrix)
Wraps the coordinates of point x
such that the returning coordinates are in the first unit cell with all-positive coordinates. The unit cell has to be a matrix of size (N,N)
.
Example
julia> using MolSimToolkitShared: wrap_to_first
julia> uc = [10.0 0.0 0.0; 0.0 10.0 0.0; 0.0 0.0 10.0]
3×3 Matrix{Float64}:
10.0 0.0 0.0
0.0 10.0 0.0
0.0 0.0 10.0
julia> wrap_to_first([15.0, 13.0, 2.0], uc)
3-element Vector{Float64}:
5.0
3.0000000000000004
2.0