Help entries
ComplexMixtures.CMTypes — TypeInternal structure or function, interface may change.
Union of types to define an isaprox for testing.
ComplexMixtures.ChemFile — Typestruct ChemFile{T<:(AbstractVector)} <: TrajectoryStructure to contain a trajectory as read by Chemfiles.jl
filename::Stringformat::AbstractStringstream::ComplexMixtures.Stream{<:Chemfiles.Trajectory}nframes::Int64sides::Vector{T} where T<:(AbstractVector)solute::Selectionsolvent::Selectionx_solute::Vector{T} where T<:(AbstractVector)x_solvent::Vector{T} where T<:(AbstractVector)natoms::Int64
ComplexMixtures.ChemFile — MethodChemFile(filename::String, solute::Selection, solvent::Selection;format="" , T::Type = SVector{3,Float64})Function open will set up the IO stream of the trajectory, fill up the number of frames field and additional parameters if required.
ComplexMixtures.Density — Typemutable struct DensityStructure to contain the density values obtained from the calculation.
solute::Float64: Default: 0.0solvent::Float64: Default: 0.0solvent_bulk::Float64: Default: 0.0
ComplexMixtures.MinimumDistance — Typestruct MinimumDistanceInternal structure or function, interface may change.
Extended help
This structure contains the information, for each molecule, of if it is within the cutoff distance of the solute, the atom indexes of the associated minimum distance, the distance, and a label to mark if the reference atom of the molecule is within the cutoff distance of the solute.
The lists of minimum-distances are stored in arrays of type Vector{MinimumDistance}. The index of this vector corresponds to the index of the molecule in the original array.
within_cutoff::Booli::Int64j::Int64d::Float64ref_atom_within_cutoff::Boold_ref_atom::Float64
ComplexMixtures.NamdDCD — Typestruct NamdDCD{T<:(AbstractVector)} <: TrajectoryStructure to contain the data of a trajectory in NAMD/DCD format.
filename::Stringstream::ComplexMixtures.Stream{<:FortranFiles.FortranFile}nframes::Int64sides::Vector{T} where T<:(AbstractVector)solute::Selectionsolvent::Selectionx_solute::Vector{T} where T<:(AbstractVector)x_solvent::Vector{T} where T<:(AbstractVector)sides_in_dcd::Boollastatom::Int64sides_read::Vector{Float64}x_read::Vector{Float32}y_read::Vector{Float32}z_read::Vector{Float32}
ComplexMixtures.NamdDCD — MethodNamdDCD(filename::String, solute::Selection, solvent::Selection;T::Type = SVector{3,Float64})This function initializes the structure above, returning the data and the vectors with appropriate lengths.
ComplexMixtures.Options — Typestruct OptionsStructure that contains the detailed input options.
firstframe::Int64: Default: 1lastframe::Int64: Default: -1stride::Int64: Default: 1irefatom::Int64: Default: -1n_random_samples::Int64: Default: 10binstep::Float64: Default: 0.02dbulk::Float64: Default: 10.0cutoff::Float64: Default: 10.0usecutoff::Bool: Default: falselcell::Int64: Default: 1GC::Bool: Default: trueGC_threshold::Float64: Default: 0.1seed::Int64: Default: 321StableRNG::Bool: Default: falsenthreads::Int64: Default: 0silent::Bool: Default: false
ComplexMixtures.OutputFiles — TypeInternal structure or function, interface may change.
mutable struct OutputFilesStructure to contain the names of the output files.
output::Stringsolute_atoms::Stringsolvent_atoms::String
ComplexMixtures.Overview — TypeInternal structure or function, interface may change.
mutable struct OverviewStructure that is used to dispatch the show of a overview.
R::Resultdomain_molar_volume::Float64: Default: 0.0density::ComplexMixtures.Density: Default: Density()solvent_molar_volume::Float64: Default: 0.0solvent_molar_volume_bulk::Float64: Default: 0.0solute_molar_volume::Float64: Default: 0.0
ComplexMixtures.PDBTraj — Typestruct PDBTraj{T<:(AbstractVector)} <: TrajectoryStructure to contain PDB trajectories. Frames must be separated by "END", and with periodic cell sizes in the "CRYST1" field.
This structure and functions can be used as a template to implement the reading of other trajectory formats.
filename::Stringstream::ComplexMixtures.Stream{<:IOStream}nframes::Int64sides::Vector{T} where T<:(AbstractVector)solute::Selectionsolvent::Selectionx_solute::Vector{T} where T<:(AbstractVector)x_solvent::Vector{T} where T<:(AbstractVector)
ComplexMixtures.PDBTraj — MethodPDBTraj(pdbfile::String, solute::Selection, solvent::Selection;T::Type = SVector{3,Float64})Function open will set up the IO stream of the trajectory, fill up the number of frames field and additional parameters if required
ComplexMixtures.Result — Typemutable struct Result{T<:VecOrMat{Float64}}Structure to contain the results of the MDDF calculation.
nbins::Int64dbulk::Float64cutoff::Float64d::Vector{Float64}: Default: zeros(nbins)md_count::Vector{Float64}: Default: zeros(nbins)md_count_random::Vector{Float64}: Default: zeros(nbins)coordination_number::Vector{Float64}: Default: zeros(nbins)coordination_number_random::Vector{Float64}: Default: zeros(nbins)mddf::Vector{Float64}: Default: zeros(nbins)kb::Vector{Float64}: Default: zeros(nbins)autocorrelation::Boolsolvent::ComplexMixtures.SolSummarysolute::ComplexMixtures.SolSummarysolute_atom::VecOrMat{Float64}: Default: zeros(nbins, solute.natomspermol)solvent_atom::VecOrMat{Float64}: Default: zeros(nbins, solvent.natomspermol)rdf_count::Vector{Float64}: Default: zeros(nbins)rdf_count_random::Vector{Float64}: Default: zeros(nbins)sum_rdf_count::Vector{Float64}: Default: zeros(nbins)sum_rdf_count_random::Vector{Float64}: Default: zeros(nbins)rdf::Vector{Float64}: Default: zeros(nbins)kb_rdf::Vector{Float64}: Default: zeros(nbins)density::ComplexMixtures.Density: Default: Density()volume::ComplexMixtures.Volume: Default: Volume(nbins)options::Optionsirefatom::Int64lastframe_read::Int64nframes_read::Int64files::Vector{String}weights::Vector{Float64}
The Result{Vector{Float64}} parametric type is necessary only for reading the JSON3 saved file.
ComplexMixtures.Selection — Typestruct SelectionStructure that contains the information about the solute and solvent molecules.
natoms::Int64nmols::Int64natomspermol::Int64index::Vector{Int64}imol::Vector{Int64}names::Vector{String}
ComplexMixtures.SolSummary — TypeInternal structure or function, interface may change.
struct SolSummaryStructures to contain the details of a solute or solvent to store in the results of the MDDF calculation.
natoms::Int64nmols::Int64natomspermol::Int64
ComplexMixtures.Trajectory — TypeTrajectory(filename::String, solute::Selection, solvent::Selection; format::String = "", chemfiles = false)Trajectory constructor data type.
Defaults to reading with the Chemfiles infrastructure, except for DCD and PDB trajectory files, if the "PDBTraj" option is provided.
See memory issue (https://github.com/chemfiles/Chemfiles.jl/issues/44)
ComplexMixtures.Units — TypeInternal structure or function, interface may change.
struct UnitsUnit conversions.
mole::Any: Default: 6.022140857e23Angs3tocm3::Any: Default: 1.0e24Angs3toL::Any: Default: 1.0e27Angs3tocm3permol::Any: Default: mole / Angs3tocm3Angs3toLpermol::Any: Default: mole / Angs3toLSitesperAngs3tomolperL::Any: Default: Angs3toL / mole
ComplexMixtures.Volume — Typemutable struct VolumeStructures to contain the volumes obtained from calculations.
total::Float64bulk::Float64domain::Float64shell::Vector{Float64}
Base.isapprox — MethodBase.isapprox(r1::T, r2::T; debug=false) where T <: CMTypesInternal structure or function, interface may change.
Function to test if two runs offered similar results. Mostly used in the package testing routines.
Base.merge — Methodmerge(r::Vector{Result})This function merges the results of MDDF calculations obtained by running the same analysis on multiple trajectories, or multiple parts of the same trajectory. It returns a Result structure of the same type, with all the functions and counters representing averages of the set provided weighted by the number of frames read in each Result set.
Base.write — Methodwrite(R::ComplexMixtures.Result, filename::String, solute::Selection, solvent::Selection)Function to write the final results to output files as simple tables that are human-readable and easy to analyze with other software
If the solute and solvent selections are provides, pass on the atom names.
Base.write — Methodwrite(R::ComplexMixtures.Result, filename::String;
solute_names::Vector{String} = ["nothing"],
solvent_names::Vector{String} = ["nothing"])Optional passing of atom names.
ComplexMixtures.VMDselect — MethodVMDselect(inputfile::String, selection::String; vmd="vmd" )Select atoms using vmd selection syntax, with vmd in background
Returns the list of index (one-based) and atom names
Function to return the selection from a input file (topology, coordinates, etc), by calling VMD in the background.
ComplexMixtures.contributions — Methodcontributions(s::Selection, atom_contributions::Matrix{Float64}, selection)Extract the contribution of a given atom type selection from the solute or solvent atomic contributions to the MDDF.
s here is the solute or solvent selection (type ComplexMixtures.Selection) atom_contributions is the R.solute_atom or R.solvent_atom arrays of the Result structure, and the last argument is the selection of atoms from the solute to be considered, given as a list of indexes, list of atom names, vector of PDBTools.Atoms, or a PDBTools.Residue.
Extended help
For selections of one molecule, the function has an additional keyword option first_atom_is_ref that is false by default. If set to true, the index first atom of the selection is considered as a reference atom. For example if a solute has 100 atoms, but its first atom in the PDB file is number 901, the selection of indexes [1, 2, 3] will refer to atoms with indexes [901, 902, 903].
ComplexMixtures.coordination_number — Functioncoordination_number(trajectory::Trajectory, options::Options)Computes the coordination numbers for each solute molecule in the trajectory, given the Trajectory. This is an auxiliary function of the ComplexMixtures package, which is used to compute coordination numbers when the normalization of the distribution is not possible or needed.
The output is a Result structure, which contains the data as the result of a call to mddf, except that all counters which require normalization of the distribution will be zero. In summary, this result data structure can be used to compute the coordination numbers, but not the MDDF, RDF, or KB integrals.
Examples
julia> trajectory = Trajectory("./trajectory.dcd",solute,solvent);
julia> results = mddf(trajectory);
julia> coordination_numbers = coordination_number(trajectory);ComplexMixtures.coordination_number — Functioncoordination_number(R::Result) = R.coordination_number
coordination_number(R::Result, group_contributions::Vector{Float64})
coordination_number(s::Selection, atom_contributions::Matrix{Float64}, R::Result, group)Computes the coordination number of a given group of atoms from the solute or solvent atomic contributions to the MDDF. If no group is defined (first call above), the coordination number of the whole solute or solvent is returned.
If the group_contributions to the mddf are computed previously with the contributions function, the result can be used to compute the coordination number by calling coordination_number(R::Result, group_contributions).
Otherwise, the coordination number can be computed directly with the second call, where:
s is the solute or solvent selection (type ComplexMixtures.Selection)
atom_contributions is the R.solute_atom or R.solvent_atom arrays of the Result structure
R is the Result structure,
and the last argument is the selection of atoms from the solute to be considered, given as a list of indexes, list of atom names, or a selection following the syntax of PDBTools, or vector of PDBTools.Atoms, or a PDBTools.Residue
Examples
In the following example we compute the coordination number of the atoms of residue 50 (of the solute) with the solvent atoms of TMAO, as a function of the distance. Finally, we show the average number of TMAO molecules within 5 Angstroms of residue 50. The findlast(<(5), R.d) part of the code below returns the index of the last element of the R.d array that is smaller than 5 Angstroms.
Precomputing the group contributions Using the contributions function
using ComplexMixtures, PDBTools
pdb = readPDB("test/data/NAMD/structure.pdb");
R = load("test/data/NAMD/protein_tmao.json");
solute = Selection(PDBTools.select(pdb, "protein"), nmols=1);
residue50 = PDBTools.select(pdb, "residue 50");
# Compute the group contributions to the MDDF
residue50_contribution = contributions(solute, R.solute_atom, residue50);
# Now compute the coordination number
residue50_coordination = coordination_number(R, residue50_contribution)
# Output the average number of TMAO molecules within 5 Angstroms of residue 50
residue50_coordination[findlast(<(5), R.d)]Without precomputing the group_contribution
using ComplexMixtures, PDBTools
pdb = readPDB("test/data/NAMD/structure.pdb");
R = load("test/data/NAMD/protein_tmao.json");
solute = Selection(PDBTools.select(pdb, "protein"), nmols=1);
residue50 = PDBTools.select(pdb, "residue 50");
# Compute the coordination number
residue50_coordination = coordination_number(solute, R.solute_atom, R, group)
# Output the average number of TMAO molecules within 5 Angstroms of residue 50
residue50_coordination[findlast(<(5), R.d)]ComplexMixtures.coordination_number_frame! — Methodcoordination_number_frame!(R::Result, system::AbstractPeriodicSystem, buff::Buffer, options::Options, RNG)Internal structure or function, interface may change.
Computes the coordination numbers for a single frame. Modifies the data in the R (type Result) structure.
ComplexMixtures.eulermat — Methodeulermat(beta, gamma, theta, deg::String)Internal structure or function, interface may change.
This routine was added because it defines the rotation in the "human" way, an is thus used to set the position of the fixed molecules. deg can only be "degree", in which case the angles with be considered in degrees. If no deg argument is provided, radians are used.
That means: beta is a counterclockwise rotation around x axis. gamma is a counterclockwise rotation around y axis. theta is a counterclockwise rotation around z axis.
ComplexMixtures.finalresults! — Methodfinalresults!(R::Result, options::Options, trajectory::Trajectory)Internal structure or function, interface may change.
Function that computes the final results of all the data computed by averaging according to the sampling of each type of data, and converts to common units.
Computes also the final distribution functions and KB integrals.
This function modified the values contained in the R data structure.
ComplexMixtures.gr — Methodgr(R::Result) = gr(R.d,R.rdf_count,R.density.solvent_bulk,R.options.binstep)If a Result structure is provided without further details, use the rdf count and the bulk solvent density.
ComplexMixtures.gr — Methodgr(r::Vector{Float64}, count::Vector{Float64}, density::Float64, binstep::Float64)Computes the radial distribution function from the count data and the density.
This is exactly a conventional g(r) if a single atom was chosen as the solute and solvent selections.
Returns both the g(r) and the kb(r)
ComplexMixtures.grid3D — Methodgrid3D(solute,solute_atoms,mddf_result,output_file; dmin=1.5, ddax=5.0, step=0.5)This function builds the grid of the 3D density function and fills an array of mutable structures of type Atom, containing the position of the atoms of grid, the closest atom to that position, and distance.
solute is a ComplexMixtuers.Selection, defining the solute. solute_atoms is the corresponding vector of PDBTools.Atoms, and mddf_result is the result of a mddf_result calculation with the correspondign solute.
dmin and dmax define the range of distance where the density grid will be built, and step defines how fine the grid must be. Be aware that fine grids involve usually a very large (hundreds of thousands points).
All parameters can be provides as keyword parameters.
Example
julia> using ComplexMixtures, PDBTools
julia> pdb = readPDB("./system.pdb");
julia> R = ComplexMixtures.load("./results.json");
julia> protein = select(pdb,"protein");
julia> solute = ComplexMixtures.Selection(protein,nmols=1);
julia> grid = ComplexMixtures.grid3D(solute=solute, solute_atoms=protein, mddf_result=R, output_file="grid.pdb");
grid will contain a vector of Atoms with the information of the MDDF at each grid point, and the same data will be written in the grid.pdb file. This PDB file can be opened in VMD, for example, and contain in the beta field the contribution of each protein residue to the MDDF at each point in space relative to the protein, and in the occupancy field the distance to the protein. Examples of how this information can be visualized are provided in the user guide of ComplexMixtures.
ComplexMixtures.itype — Methoditype(iatom::Int, natomspermol::Int)Internal structure or function, interface may change.
Given the index of the atom in the vector of coordinates of the solute or the solvent, returns the type of the atom, that is, the index of this atom within the molecule (goes from 1 to natomspermol)
ComplexMixtures.load — Methodload(filename::String)Function to load the json saved results file into the Result data structure.
ComplexMixtures.mddf — Functionmddf(trajectory::Trajectory, options::Options)Function that computes the minimum-distance distribution function, atomic contributions, and KB integrals, given the Trajectory structure of the simulation and, optionally, parameters given as a second argument of the Options type. This is the main function of the ComplexMixtures package.
Examples
julia> trajectory = Trajectory("./trajectory.dcd",solute,solvent);
julia> results = mddf(trajectory);or, to set some custom optional parameter,
julia> options = Options(lastframe=1000);
julia> results = mddf(trajectory,options);ComplexMixtures.mddf_frame! — Methodmddf_frame!(R::Result, system::AbstractPeriodicSystem, buff::Buffer, options::Options, RNG)Internal structure or function, interface may change.
Computes the MDDF for a single frame. Modifies the data in the R (type Result) structure.
ComplexMixtures.minimum_distances! — Methodminimum_distances!(system::CellListMap.PeriodicSystem, R::Result)Internal structure or function, interface may change.
Function that computes the list of distances of solvent molecules to a solute molecule. It updates the lists of minimum distances.
ComplexMixtures.mol_index — Methodmol_index(i_atom, natomspermol) = (i_atom-1) ÷ natomspermol + 1Internal structure or function, interface may change.
Extended help
Sets the index of the molecule of an atom in the simples situation, in which all molecules have the same number of atoms.
ComplexMixtures.mol_range — Methodmol_range(imol, n_atoms_per_molecule)Internal structure or function, interface may change.
Given the index and the number of atoms per molecule, returns the range of indices of of an array of coordinates that corresponds to the molecule.
ComplexMixtures.move! — Methodmove!(x::AbstractVector, newcm::AbstractVector,beta, gamma, theta)Internal structure or function, interface may change.
Translates and rotates a molecule according to the desired input center of coordinates and Euler rotations modifyies the vector x.
ComplexMixtures.overview — Methodoverview(R::Result)Function that outputs the volumes and densities in the most natural units.
ComplexMixtures.random_move! — Methodrandom_move!(x_ref::AbstractVector{T},
irefatom::Int,
system::AbstractPeriodicSystem,
x_new::AbstractVector{T}, RNG) where {T<:SVector}Internal structure or function, interface may change.
Function that generates a new random position for a molecule.
The new position is returned in x_new, a previously allocated array.
ComplexMixtures.randomize_solvent! — Methodrandomize_solvent!(system, buff, n_solvent_in_bulk, options, RNG)Internal structure or function, interface may change.
Generate a random solvent distribution from the bulk molecules of a solvent
ComplexMixtures.save — Methodsave(R::Result, filename::String)Function to write the result data structure to a json file.
ComplexMixtures.setbin — Methodsetbin(d,step)Internal structure or function, interface may change.
Function that sets to which histogram bin a data point pertains simple, but important to keep consistency over all calls.
ComplexMixtures.setup_PeriodicSystem — Methodsetup_PeriodicSystem(trajectory::Trajectory, options::Options)Internal structure or function, interface may change.
Setup the periodic system from CellListMap, to compute minimimum distances. The system will be setup such that xpositions corresponds to one molecule of the solute, and ypositions contains all coordinates of all atoms of the solvent.
ComplexMixtures.shellradius — Methodshellradius(i,step)Internal structure or function, interface may change.
Compute the point in which the radius comprises half of the volume of the shell.
ComplexMixtures.sphereradiusfromshellvolume — Methodsphereradiusfromshellvolume(volume,step)Internal structure or function, interface may change.
Computes the radius that corresponds to a spherical shell of a given volume.
ComplexMixtures.sphericalshellvolume — Methodsphericalshellvolume(i,step)Internal structure or function, interface may change.
Computes the volume of the spherical shell defined within [(i-1)step,istep].
ComplexMixtures.sum! — Methodsum!(R1::Result, R2::Result)Internal structure or function, interface may change.
Sum the counts of two Results structures, adding the result to the first structure as in R1 = R1 + R2.
ComplexMixtures.title — Methodtitle(R::Result, solute::Selection, solvent::Selection)
title(R::Result, solute::Selection, solvent::Selection, nspawn::Int)Internal structure or function, interface may change.
Print some information about the run.
ComplexMixtures.update_list! — Methodupdate_list!(i, j, d2, iref_atom::Int, mol_index_i::F, isolute::Int, list::Vector{MinimumDistance{T}}) where {F<:Function, T}Internal structure or function, interface may change.
Function that updates a list of minimum distances given the indexes of the atoms involved for one pair within cutoff, for autocorrelations (such that the identity of isolute is needed)
ComplexMixtures.update_list! — Methodupdate_list!(i, j, d2, iref_atom::Int, mol_index_i::F, list::Vector{MinimumDistance{T}}) where {F<:Function, T}Internal structure or function, interface may change.
Function that updates a list of minimum distances given the indexes of the atoms involved for one pair within cutoff.
ComplexMixtures.update_md — Methodupdate_md(md1::MinimumDistance{T}, md2::MinimumDistance{T}) where {T}Internal structure or function, interface may change.
Function that returns the updated minimum distance structure after comparing two structures associated with the same molecule.
ComplexMixtures.updatecounters! — Methodupdatecounters!(R::Result, system::AbstractPeriodicSystem)Internal structure or function, interface may change.
Function that updates the minimum-distance counters in R
ComplexMixtures.viewmol — Methodviewmol(i::Int, x::Vector{T}, n::Int) where TInternal structure or function, interface may change.
Returns a view of a coordinate vector corresponding to the atoms of a molecule with index i. n is the number of atoms of the molecule.
ComplexMixtures.which_types — Methodwhich_types(s::Selection, indexes::Vector{Int})Internal structure or function, interface may change.
Function that returns the list of the indexes of the types of the atoms in a selection. For example, if a selection corresponds to a solvent of water molecules: There are three types, 1, 2, and 3, corresponding to the three atoms of the water molecule. If the indexes provided are, for instance, 11, 12, and 13, corresponding to a water molecule, this function will return 1, 2 and 3.
This is used to get equivalent-atom contributions to the distribution functions. For example, the input indexes span all water molecules, the output of this function will be still the three indexes corresponding to the three types of atoms that exist in a water molecule.
It is not possible to compute the contribution of one individual water molecule if the distribution function was computed for all molecules. Thus, the necessity to identify the types of atoms involved in a selection.
ComplexMixtures.writexyz — Methodwritexyz(x::Vector{T}, file::String) where T <: AbstractVectorInternal structure or function, interface may change.
Print test xyz file.