Structural alignment
Some functions to compute the RMSD and perform rigid-body structural alignments:
MolSimToolkit.rmsd_matrix
— Methodrmsd_matrix(
simulation::Simulation,
indices::AbstractVector{Int};
mass::Union{AbstractVector{Int}, Nothing} = nothing,
align::Bool = true,
show_progress = true,
)
Computes the RMSD matrix for a set of atoms along a trajectory.
The indices
vector contains the indices of the atoms to be considered. The mass
argument can be used to provide the mass of the atoms if they are not the same. The align
argument can be used to align the frames before computing the RMSD.
The show_progress
argument can be used to show a progress bar.
Returns
A symetric matrix with the RMSD values between each pair of frames. For example, in a trajectory with 5 frames, the matrix will be a 5x5 matrix with the RMSD values between the structures of each pair of frames.
Example
julia> using MolSimToolkit, MolSimToolkit.Testing
julia> using PDBTools
julia> atoms = readPDB(Testing.namd_pdb);
julia> cas = findall(Select("name CA"), atoms); # CA indices
julia> simulation = Simulation(Testing.namd_pdb, Testing.namd_traj);
julia> rmsd_matrix(simulation, cas; show_progress=false)
5×5 Matrix{Float64}:
0.0 2.83887 2.9777 2.46214 3.80357
2.83887 0.0 2.35492 2.64463 4.68028
2.9777 2.35492 0.0 2.08246 3.46149
2.46214 2.64463 2.08246 0.0 2.97835
3.80357 4.68028 3.46149 2.97835 0.0
MolSimToolkitShared.rmsd
— Methodrmsd(simulation::Simulation, indices::AbstractVector{Int}; mass = nothing, reference_frame = nothing, show_progress = true)
Computes the root mean square deviation (RMSD) between two sets of points in along a trajectory.
Arguments
indices
vector contains the indices of the atoms to be considered.mass
argument can be used to provide the mass of the atoms if they are not the same.reference_frame
argument can be used to provide a reference frame to align the trajectory to:- If
reference_frame == nothing
, the first frame will be used (default behavior). - If
reference_frame == :average
, the average structure will be used. - If
reference_frame
is an integer, the frame at that index will be used as reference.
- If
Examples
Computing the rmsd along a trajectory
julia> using MolSimToolkit, MolSimToolkit.Testing
julia> using PDBTools
julia> atoms = readPDB(Testing.namd_pdb);
julia> simulation = Simulation(Testing.namd_pdb, Testing.namd_traj);
julia> cas = findall(sel"name CA", atoms); # CA indices
julia> rmsd(simulation, cas; show_progress=false)
5-element Vector{Float64}:
0.0
2.8388710154609034
2.9776998440690385
2.4621444212469483
3.8035683196100796
julia> rmsd(simulation, cas; reference_frame=:average, show_progress=false)
5-element Vector{Float64}:
1.8995986972454748
2.1512244220536973
1.5081703191869376
1.1651111324544219
2.757039151265317
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