Functions

Index

File IO

MolHandler.read_dcdFunction
read_dcd(filename::String;
         frame_indices::Union{Vector, OrdinalRange, Colon, Nothing} = nothing,
         step::Union{Int, Nothing} = nothing)
::Trajectory

Return Trajectory object which filled coordinates, nframe, natom fields. If you set frame_indices, only specified frame is read. If you set step, every step frame is read. You can't specify both frame_indices and step at the same time.

source
MolHandler.write_dcdFunction
write_dcd(filename::String, trj::Trajectory;
          save_step::Integer = 1, total_step::Integer = trj.nframe,
          unit_num::Integer  = 1,  time_step::Real = 1.0f0)

Write coordinates, nframe and natom information of trj to dcd file, named filename. dcd file contain other information - save interval step, unit(chain) number in the system, total step and time step of original trajectory -, and if you do not specify these parameter, these are filled with default value.

source
MolHandler.read_pdbFunction
read_pdb(filename::AbstractString; model = :unspecified)::Trajectory

Return Trajectory object which filled all field. If you set :AA to model field, this function set the mass of particle based on atomname field. If you set :CA to model field, this function only read CA atom and set the mass of particcle mase on resname field.

source
MolHandler.write_pdbFunction
write_pdb(filename::AbstractString, trj::Trajectory;
          tempfactor::Union{RealT, Nothing} = nothing) where RealT <: Real

Write coordinates, resname, resid, atomname, atomid to filename based on trj. If you set Real value to tempfactor, all tempfactor will be field with that value.

source
MolHandler.read_xyzFunction
read_xyz(filename::AbstractString)::Trajectory

Return Trajectory object which filled coordinates, nframe, natom fields.

source
MolHandler.write_xyzFunction
write_xyz(filename:AbstractString, trj::Trajectory)

Write coordinates, nframe and natom information of trj to xyz file, named filename.

source

Trajectory handling

MolHandler.get_frameFunction
get_frame(frame_idx::Int64, trajectory::Trajectory)
::Frame

Return Frame object which correspond to framd_idx frame from trajectory.

source
MolHandler.get_atomFunction
get_atom(query::Int64, trajectory::Trajectory)
::Vector{Atom}

Return Vector of queryth Atom object from trajectory.

source
get_atom(query::Int64, frame::Frame)
::Atom

Return queryth Atom object from frame.

source
MolHandler.clip_trajectoryFunction
clip_trajectory(query::Union{Integer, Vector{IntT}, OrdinalRange}, trj::Trajectory;
                query_key::Symbol = :frame)
::Trajecotory where IntT <: Integer

Clip a part of trajectory you specified. Expected query key is one of the following

  • :frame : in this case, query specify the frame index.
  • :atom : in this case, query specify the atom index.
source
MolHandler.geometric_center_of_massFunction
geometric_center_of_mass(coordinates::Vector{Coordinate};
                         atom_indices::Union{Vector, OrdinalRange, Colon} = :)
::Coordinate

Calculate the geometric center of mass of Coordinate vector for specified atom indices.

For trajectory case, geometric center of mass is calculated by center_of_mass method with geometric flag.

source
MolHandler.center_of_massFunction
center_of_mass(coordinates::Vector{Coordiante}, mass_vec{Real};
               atom_indices::Union{Vector, OrdinalRange, Colon} = :)
::Coordinate{RealT}

Calculate the (non-geometric) center of mass of Coordinate vector for specified atom indices. You have to pass mass information vector which have same length to coordinates vector.

source
center_of_mass(query::Trajectory;
               frame_indices::Union{Vector, OrdinalRange, Colon} = :,
               atom_indices::Union{Vector, OrdinalRange, Colon} = :,
               geometric::Bool = false)
::Vector{Coordinate{Real}}

Calculate the center of mass of trajectory for specified frame and atom indices. If you set geometric is true, this function calculate geometric center of mass.

source
MolHandler.move_pbc_centerFunction
move_pbc_center(coordinates::Vector{Coordinate},
                new_center::Coordinate, box_size::Coordinate)
                )::Vector{Coordinate}

Fix coordinate to consistent with new box center and box size. Box size should match to original periodic boundary box size.

source
MolHandler.radius_of_gyrationFunction
radius_of_gyration(trj::Trajectory;
                   frame_indices::Union{Vector, OrdinalRange, Colon} = :,
                   atom_indices ::Union{Vector, OrdinalRange, Colon} = :,
                   geometric::Bool = false)
::Vector{Real}

Calculate radius of gyration over the trajectory. If you set geometric is true, this function calculate radius of gyration without mass info.

source
MolHandler.pair_length_matrixFunction
pair_length_matrix(corrds1::Vector{Coordinate}, coords2::Vector{Coordinate})
::Matrix{Real}

Calculate distance matrix for all combination between coords1 and coord2. The row of returned matrix coresspond to coords1, and the column correspond to coords2.

source
pair_length_matrix(trj::Trajectory;
                   frame_indices::Union{Vector, OrdinalRange, Colon}       = :,
                   first_atom_indices::Union{Vector, OrdinalRange, Colon}  = :,
                   second_atom_indices::Union{Vector, OrdinalRange, Colon} = atom_indices1)
::Vector{Matrix{Coordinate}}

Calculate distance matrix for all combination between first_atom_indices and second_atom_indices. This cauculation apply to each frame of trajectory and the result matrices are stored in Vector. The target frame can be restricted by pass indeces vector or range to frame_indices.

source
MolHandler.pair_length_matrix_pbcFunction
 pair_length_matrix_pbc(corrds1::Vector{Coordinate}, coords2::Vector{Coordinate},
                        box_size::Coordinate)
::Matrix{Real}

Calculate distance matrix for all combination between coords1 and coord2 considering periodic boundary condition. The box size information used for calculation is specified by box_size argument. The row of returned matrix coresspond to coords1, and the column correspond to coords2.

source
pair_length_matrix_pbc(trj::Trajectory, box_size::Coordinate;
                       frame_indices::Union{Vector, OrdinalRange, Colon}       = :,
                       first_atom_indices::Union{Vector, OrdinalRange, Colon}  = :,
                       second_atom_indices::Union{Vector, OrdinalRange, Colon} = atom_indices1)
::Vector{Matrix{Coordinate}}

Calculate distance matrix for all combinations between first_atom_indices and second_atom_indices considering periodic boundary condition. The box size information used for calculation is specified by box_size argument. This cauculation apply to each frame of trajectory and the result matrices are stored in Vector. The target frame can be restricted by pass indeces vector or range to frame_indices.

source
MolHandler.contact_bool_matrixFunction
contact_bool_matrix(threshold::Real, trj::Trajectory;
                    frame_indices::Union{Array, OrdinalRange, Colon}       = :,
                    first_atom_indices::Union{Array, OrdinalRange, Colon}  = :,
                    second_atom_indices::Union{Array, OrdinalRange, Colon} = first_atom_indices)
::Vector{Matrix{Bool}}

Judge contact is formed or not. If the distance between two coordinate is shorter than threshold, contact is considered to be formed. In returned vector of matrices, each matrix correspond to contact matrix of each frame. You can specify the target frames or atoms by frame_indices, first_atom_indices or second_atom_indices. When you specify the target atoms, the row of matrices corresponds to firstatomindices and column of matrices corresponds to secondatomindices.

source
MolHandler.contact_bool_matrix_pbcFunction
contact_bool_matrix_pbc(threshold::Real, trj::Trajectory, box_size::Coordinate;
                        frame_indices::Union{Array, OrdinalRange, Colon}       = :,
                        first_atom_indices::Union{Array, OrdinalRange, Colon}  = :,
                        second_atom_indices::Union{Array, OrdinalRange, Colon} = first_atom_indices)
::Vector{Matrix{Bool}}

Judge contact is formed or not considering periodic boundary condition. If the distance between two coordinate is shorter than threshold, contact is considered to be formed. In returned vector of matrices, each matrix correspond to contact matrix of each frame. The box size information used for calculation is specified box_size argument. You can specify the target frames or atoms by frame_indices, first_atom_indices or second_atom_indices. When you specify the target atoms, the row of matrices corresponds to firstatomindices and column of matrices corresponds to secondatomindices.

source
MolHandler.contact_probability_matrixFunction
contact_probability_matrix(threshold::Real, trj::Trajectory;
                           frame_indices::Union{Vector, OrdinalRange, Colon}       = :,
                           first_atom_indices::Union{Vector, OrdinalRange, Colon}  = :,
                           second_atom_indices::Union{Vector, OrdinalRange, Colon} = first_atom_indices)
::Matrix{Real}

Calculate contact formation probability over the trajectory. If the distance between two coordinate is shorter than threshold, contact is considered to be formed. You can specify the target frames or atoms by frame_indices, first_atom_indices or second_atom_indices. When you specify the target atoms, the row of matrices corresponds to firstatomindices and column of matrices corresponds to secondatomindices.

source
MolHandler.fix_pbcFunction
fix_pbc(coordinates::Vector{Coordinate}, groupid_vec::Vector{Integer},
box_size::Coordinate{<:Real};
x::Bool = true, y::Bool = true, z::Bool = true)::Vector{Coordinate}

Fix the atom group splited by periodic boundary box. This group is defined by groupid_vec which contain each particles group id. For example, if your system have 3 atom, and atom 1 and 2 are group 1, and atom 3 is group 2, this groupid_vec is [1, 1, 2]. If the group over the boundary of the box, the atoms in the group which separated from first atom of that move to the position where the atom locate without periodic boundary condition. You can specify the axies which applied fixing by x,y and option. If you set z = false, only x and y coordinate fixed and z don't change.

source
fix_pbc(trj::Trajectory, groupid_vec::Vector{Integer},
box_size::Coordinate;
x::Bool = true, y::Bool = true, z::Bool = true)::Trajectory

Fix residues splited by periodic boundary condition. This is more specific version of fix_pbc function for trajectory handling. For example, if your system have 3 atom, and atom 1 and 2 are group 1, and atom 3 is group 2, this groupid_vec is [1, 1, 2]. If the group over the boundary of the box, the atoms in the group which more than half a box away from the first atom of that move to the position where the atom locate without periodic boundary condition. You can specify the axies which applied fixing by x,y and option. If you set z = false, only x and y coordinate fixed and z don't change.

source

Structure analysis

MolHandler.sasaFunction
sasa(coordinates::Vector{Coordinate}, radiuses::Vector{Real},
     solvent_radius::RealT = 1.4,     MCtrial_num::Integer = 1000
     seed::Integer         = undef)::Vector{Real}

Calculate Solvent Accessible Surface Area(SASA) for each particle. The length of radiuses should match to the number of particles in coordinates. The default value of solvent radius, 1.4, corresponds to water case 1.4 Å. MCtrial means MonteCarlo trial num for calculate the solvent accessible area of each particle. For each particle, this function generate the point MCtrial times on the solvent accessible surface and count the number of non-overlapping dots with other particle.

source

Utility

MolHandler.atom_massFunction
atom_mass(name::AbstractString)
::Float32

Return mass of the name atom. The kind of each atom is judged by PDB format atom names rule. https://cdn.rcsb.org/wwpdb/docs/documentation/file-format/PDBformat1992.pdf

source
MolHandler.residue_massFunction
residue_mass(name::AbstractString)
::Float32

Return mass of the name residue. The kind of each residue is judged by amino acid 3-letter abbreviation.

source
MolHandler.distance_pbcFunction
distance_pdc(first_atom::Coordinate,  second_atom::Coordinate,
             box_size::Coordinate)
::Real

Calculate distance first_atom and second_atom considering periodic boundary condition. The box size information used calculatin is specified by box_size argument.

source

Multi-Threading

MolHandler.pair_length_matrix_parallelFunction
pair_length_matrix_parallel(trj::Trajectory;
                            frame_indices::Union{Vector, OrdinalRange, Colon}       = :,
                            first_atom_indices::Union{Vector, OrdinalRange, Colon}  = :,
                            second_atom_indices::Union{Vector, OrdinalRange, Colon} = atom_indices1)
::Vector{Matrix{Coordinate}}

Multi-threads version of pairlengthmatrix. If you set available threads number to Threads.nthreads(), this function would faster than non-parallel version.

source
MolHandler.contact_bool_matrix_parallelFunction
contact_bool_matrix_parallel(threshold::Real, trj::Trajectory;
                             frame_indices::Union{Vector, OrdinalRange, Colon}       = :,
                             first_atom_indices::Union{Vector, OrdinalRange, Colon}  = :,
                             second_atom_indices::Union{Vector, OrdinalRange, Colon} = first_atom_indices)
::Vector{Matrix{Bool}}

Multi-threads version of contactboolmatrix. If you set available threads number to Threads.nthreads(), this function would faster than non-parallel version.

source
MolHandler.contact_probability_matrix_parallelFunction
contact_probability_matrix_parallel(threshold::Real, trj::Trajectory;
                                    frame_indices::Union{Vector, OrdinalRange, Colon}       = :,
                                    first_atom_indices::Union{Vector, OrdinalRange, Colon}  = :,
                                    second_atom_indices::Union{Vector, OrdinalRange, Colon} = first_atom_indices)
::Matrix{Real}

Multi-threads version of contactprobabilitymatrix. If you set available threads number to Threads.nthreads(), this function would faster than non-parallel version.

source