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
positions
coordination_number, bulk_coordination
center_of_mass
wrap, wrap_to_first
align, align!, rmsd, alignment_movements, apply_alignment_transformations!
dihedral, dihedrals
get_atomsMethods
MolSimToolkitShared.align — Function
align(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! — Function
align(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.alignment_movements — Method
alignment_movements(
x::AbstractVector{<:AbstractVector},
y::AbstractVector{<:AbstractVector};
mass=nothing,
xm=zeros(3, length(x)),
xp=zeros(3, length(x)),
)Returns cmx, cmy, and rotation_matrix, necessary to translate and rotate y to minimize the RMSD relative to x.
Can be used to obtain the movement that aligns two sets of points, of the same size, and apply these movements to another set of points.
MolSimToolkitShared.apply_alignment_transformation! — Method
apply_alignment_transformation!(x, cmx, cmy, rotation_matrix)Translates and rotates x such that the coordinates cmx are moved to cmy, and such apply the rotation defined in rotation_matrix.
MolSimToolkitShared.center_of_mass — Method
center_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.dihedral — Method
dihedral(v1, v2, v3, v4; degrees=true)
dihedral(v::AbstractVector; degrees=true)Computes the dihedral angle between the planes formed by the vectors v1-v2 and v2-v3, and v2-v3 and v3-v4. The input vectors must have 3 elements. The function returns the dihedral angle in radians or degrees.
If the input is a vector with 4 vectors, the function computes the dihedral angle between the planes formed by the vectors v[1]-v[2] and v[2]-v[3], and v[2]-v[3] and v[3]-v[4].
The optional argument degrees specifies whether the output is in degrees (default) or radians.
MolSimToolkitShared.dihedrals — Method
dihedrals(v::AbstractVector{<:AbstractVector}; degrees=true)Computes the dihedral angles for many sets of 4 vectors. The input is a vector of vectors, where each element is a vector with 4 vectors. The function returns a vector with the dihedral angles in radians or degrees.
Example
julia> using MolSimToolkitShared: dihedrals
julia> v1 = [[-8.483, -14.912, -6.726], [-5.113, -13.737, -5.466], [-3.903, -11.262, -8.062], [-1.162, -9.64, -6.015]];
julia> v2 = [[-9.229, -14.861, -5.481], [-8.483, -14.912, -6.726], [-7.227, -14.047, -6.599], [-7.083, -13.048, -7.303]];
julia> dihedrals([v1,v2])
2-element Vector{Float64}:
164.43481280739516
-115.82544005374316MolSimToolkitShared.rmsd — Method
rmsd(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.0MolSimToolkitShared.wrap — Function
wrap(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 — Function
wrap_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