IMP
2.3.0
The Integrative Modeling Platform
|
Functionality for loading, creating, manipulating and scoring atomic structures. More...
Functionality for loading, creating, manipulating and scoring atomic structures.
Molecules and collections of molecules are represented using Hierarchy particles. Whenever possible, one should prefer to use IMP::atom::Selection and related helper functions to manipulate molecules as that provides a simple and biologically relevant way of describing parts of molecules of interest.
The name "residue index" is used to refer to the index of the residue in the conventional description of the protein, as opposed to its index among the set of residues which are found in the current molecule. The same concept is known as the "residue number" in PDB files. This index is not necessarily unique within a Chain; however, the combination of the residue index and insertion code should be.
Author(s): Daniel Russel, Ben Webb, Dina Schneidman, Javier Velazquez-Muriel
Maintainer: benmwebb
License: LGPL This library is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version.
Publications:
Classes | |
class | AllMol2Selector |
Read all atoms. More... | |
class | AllPDBSelector |
Defines a selector that will pick every ATOM and HETATM record. More... | |
class | AndPDBSelector |
Select atoms which are selected by both selectors. More... | |
class | Angle |
A particle that describes an angle between three particles. More... | |
class | AngleSingletonScore |
Score the angle based on a UnaryFunction,. More... | |
class | Atom |
A decorator for a particle representing an atom. More... | |
class | ATOMPDBSelector |
Select all non-alternative ATOM records. More... | |
class | AtomType |
The type of an atom. More... | |
class | BackbonePDBSelector |
Select all backbone (N,CA,C,O) ATOM records. More... | |
class | BerendsenThermostatOptimizerState |
Maintains temperature during molecular dynamics. More... | |
class | Bond |
A decorator for wrapping a particle representing a molecular bond. More... | |
class | Bonded |
A decorator for a particle which has bonds. More... | |
class | BondedPairFilter |
A filter for bonds. More... | |
class | BondEndpointsRefiner |
Return the endpoints of a bond. More... | |
class | BondGeometry |
Display a Bond particle as a segment. More... | |
class | BondGraph |
Represent a bond graph as a boost graph. More... | |
class | BondPairContainer |
A container that returns pairs of the endpoints of the bonds. More... | |
class | BondsGeometry |
Display an IMP::SingletonContainer of Bond particles as segments. More... | |
class | BondSingletonScore |
Score the bond based on a UnaryFunction,. More... | |
class | BrownianDynamics |
Simple Brownian dynamics simulator. More... | |
class | CAlphaPDBSelector |
Select all CA ATOM records. More... | |
class | CBetaPDBSelector |
Select all CB ATOM records. More... | |
class | CenterOfMass |
A particle that is the center of mass of other particles. More... | |
class | Chain |
Store info for a chain of a protein. More... | |
class | ChainPDBSelector |
Select all ATOM and HETATM records with the given chain ids. More... | |
class | Charged |
A decorator for a point particle that has an electrostatic charge. More... | |
class | CHARMMAtom |
A decorator for an atom that has a defined CHARMM type. More... | |
class | CHARMMAtomTopology |
A single atom in a CHARMM topology. More... | |
class | CHARMMBondEndpoint |
The end of a bond, angle, dihedral, improper, or internal coordinate. More... | |
struct | CHARMMBondParameters |
The parameters for a CHARMM bond or angle. More... | |
class | CHARMMConnection |
A connection (bond, angle, dihedral) between some number of endpoints. More... | |
struct | CHARMMDihedralParameters |
The parameters for a CHARMM dihedral or improper. More... | |
class | CHARMMIdealResidueTopology |
The ideal topology of a single residue. More... | |
class | CHARMMInternalCoordinate |
A geometric relationship between four atoms. More... | |
class | CHARMMParameters |
CHARMM force field parameters. More... | |
class | CHARMMPatch |
A CHARMM patch residue. More... | |
class | CHARMMResidueTopology |
The topology of a single residue in a model. More... | |
class | CHARMMResidueTopologyBase |
Base class for all CHARMM residue-based topology. More... | |
class | CHARMMSegmentTopology |
The topology of a single CHARMM segment in a model. More... | |
class | CHARMMStereochemistryRestraint |
Enforce CHARMM stereochemistry on the given Hierarchy. More... | |
class | CHARMMTopology |
The topology of a complete CHARMM model. More... | |
class | Copy |
A decorator for keeping track of copies of a molecule. More... | |
class | CoulombPairScore |
Coulomb (electrostatic) score between a pair of particles. More... | |
class | CoverBond |
Cover a bond with a sphere. More... | |
class | CPDBSelector |
Select all C (not CA or CB) ATOM records. More... | |
class | Diffusion |
A decorator for a diffusing particle with a diffusion coefficient. More... | |
class | Dihedral |
A particle that describes a dihedral angle between four particles. More... | |
class | DihedralSingletonScore |
Score the dihedral angle. More... | |
class | Domain |
A decorator to associate a particle with a part of a protein. More... | |
class | DopePairScore |
class | ElementTable |
class | EzRestraint |
Ez Potential kernel::Restraint. More... | |
class | ForceFieldParameters |
Storage and access to force field. More... | |
class | ForceSwitch |
Smooth interaction scores by switching the derivatives (force switch). More... | |
class | Fragment |
A decorator to associate a particle with a part of a protein/DNA/RNA. More... | |
class | HierarchiesGeometry |
Display an IMP::SingletonContainer of IMP::atom::Hierarchy particles as balls. More... | |
class | Hierarchy |
The standard decorator for manipulating molecular structures. More... | |
class | HierarchyGeometry |
Display an IMP::atom::Hierarchy particle as balls. More... | |
class | HydrogenPDBSelector |
Select all hydrogen ATOM and HETATM records. More... | |
class | ImproperSingletonScore |
Score the improper dihedral based on a UnaryFunction,. More... | |
class | LangevinThermostatOptimizerState |
Maintains temperature during molecular dynamics. More... | |
class | LennardJones |
A decorator for a particle that has a Lennard-Jones potential well. More... | |
class | LennardJonesPairScore |
Lennard-Jones score between a pair of particles. More... | |
class | Mass |
Add mass to a particle. More... | |
class | Mol2Selector |
A base class for choosing which Mol2 atoms to read. More... | |
class | MolecularDynamics |
Simple molecular dynamics optimizer. More... | |
class | Molecule |
A decorator for a molecule. More... | |
class | NonAlternativePDBSelector |
Select all ATOM and HETATM records which are not alternatives. More... | |
class | NonHydrogenMol2Selector |
Defines a selector that will pick only non-hydrogen atoms. More... | |
class | NonWaterNonHydrogenPDBSelector |
Select non water and non hydrogen atoms. More... | |
class | NonWaterPDBSelector |
Select all non-water non-alternative ATOM and HETATM records. More... | |
class | NotPDBSelector |
Select atoms which are not selected by a given selector. More... | |
class | NPDBSelector |
Select all N ATOM records. More... | |
class | OrientedSoapPairScore |
Score a pair of atoms using an orientation-dependent SOAP score. More... | |
class | OrPDBSelector |
Select atoms which are selected by either selector. More... | |
class | PDBSelector |
Select which atoms to read from a PDB file. More... | |
class | PPDBSelector |
Select all P (= phosphate) ATOM records. More... | |
class | ProteinLigandAtomPairScore |
class | ProteinLigandRestraint |
Score a pair of molecules. More... | |
class | RemoveRigidMotionOptimizerState |
Removes rigid translation and rotation from the particles. More... | |
class | RemoveTranslationOptimizerState |
Removes rigid translation from the particles. More... | |
class | Representation |
A decorator for a representation. More... | |
class | Residue |
A decorator for a residue. More... | |
class | ResidueType |
The type for a residue. More... | |
class | RigidBodyDiffusion |
class | SameResiduePairFilter |
class | SecondaryStructureResidue |
A decorator for a residue with probability of secondary structure. More... | |
class | Selection |
Select hierarchy particles identified by the biological name. More... | |
class | SelectionGeometry |
Display a Selection. More... | |
class | Simulator |
The base class for simulators. More... | |
class | SmoothingFunction |
Base class for smoothing nonbonded interactions as a function of distance. More... | |
class | SoapPairFilter |
Filter atom pairs for SOAP. More... | |
class | State |
Associate an integer "state" index with a hierarchy node. More... | |
class | StereochemistryPairFilter |
A filter that excludes bonds, angles and dihedrals. More... | |
class | VelocityScalingOptimizerState |
Maintains temperature during molecular dynamics by velocity scaling. More... | |
class | WaterPDBSelector |
Select all non-water ATOM and HETATM records. More... | |
class | WritePDBOptimizerState |
Enumerations | |
enum | Element { UNKNOWN_ELEMENT = 0, H = 1, He = 2, Li = 3, Be = 4, B = 5, C = 6, N = 7, O = 8, F = 9, Ne = 10, Na = 11, Mg = 12, Al = 13, Si = 14, P = 15, S = 16, Cl = 17, Ar = 18, K = 19, Ca = 20, Sc = 21, Ti = 22, V = 23, Cr = 24, Mn = 25, Fe = 26, Co = 27, Ni = 28, Cu = 29, Zn = 30, Ga = 31, Ge = 32, As = 33, Se = 34, Br = 35, Kr = 36, Rb = 37, Sr = 38, Y = 39, Zr = 40, Nb = 41, Mo = 42, Tc = 43, Ru = 44, Rh = 45, Pd = 46, Ag = 47, Cd = 48, In = 49, Sn = 50, Sb = 51, Te = 52, I = 53, Xe = 54, Cs = 55, Ba = 56, La = 57, Ce = 58, Pr = 59, Nd = 60, Pm = 61, Sm = 62, Eu = 63, Gd = 64, Tb = 65, Dy = 66, Ho = 67, Er = 68, Tm = 69, Yb = 70, Lu = 71, Hf = 72, Ta = 73, W = 74, Re = 75, Os = 76, Ir = 77, Pt = 78, Au = 79, Hg = 80, Tl = 81, Pb = 82, Bi = 83, Po = 84, At = 85, Rn = 86, Fr = 87, Ra = 88, Ac = 89, Th = 90, Pa = 91, U = 92, Np = 93, Pu = 94, Am = 95, Cm = 96, Bk = 97, Cf = 98, Es = 99, Fm = 100, Md = 101, No = 102, Lr = 103, Db = 104, Jl = 105, Rf = 106 } |
The various elements currently supported/known. More... | |
enum | GetByType { ATOM_TYPE, RESIDUE_TYPE, CHAIN_TYPE, MOLECULE_TYPE, DOMAIN_TYPE, FRAGMENT_TYPE, XYZ_TYPE, XYZR_TYPE, MASS_TYPE, STATE_TYPE } |
enum | RepresentationType { BALLS = 0, GAUSSIANS = 1 } |
Functions | |
AtomType | add_atom_type (std::string name, Element e) |
Create a new AtomType. More... | |
void | add_bonds (Hierarchy d, const ForceFieldParameters *ffp=get_all_atom_CHARMM_parameters()) |
void | add_radii (Hierarchy d, const ForceFieldParameters *ffp=get_all_atom_CHARMM_parameters(), FloatKey radius_key=FloatKey("radius")) |
Bond | create_bond (Bonded a, Bonded b, Int t) |
Connect the two wrapped particles by a bond. More... | |
Bond | create_bond (Bonded a, Bonded b, Bond o) |
Connect the two wrapped particles by a custom bond. More... | |
Hierarchy | create_clone (Hierarchy d) |
Clone the Hierarchy. More... | |
Hierarchy | create_clone_one (Hierarchy d) |
Clone the node in the Hierarchy. More... | |
IMP::core::RigidBody | create_compatible_rigid_body (Hierarchy h, Hierarchy reference) |
Rigidify a molecule or collection of molecules. More... | |
kernel::Restraint * | create_connectivity_restraint (const Selections &s, double k, std::string name="Connectivity%1%") |
Create a restraint connecting the selections. More... | |
kernel::Restraint * | create_connectivity_restraint (const Selections &s, double x0, double k, std::string name="Connectivity%1%") |
Create a restraint connecting the selections. More... | |
core::XYZR | create_cover (const Selection &s, std::string name=std::string()) |
Bond | create_custom_bond (Bonded a, Bonded b, Float length, Float stiffness=-1) |
Connect the two wrapped particles by a custom bond. More... | |
kernel::Restraint * | create_distance_restraint (const Selection &n0, const Selection &n1, double x0, double k, std::string name="Distance%1%") |
kernel::Restraint * | create_excluded_volume_restraint (const Hierarchies &hs, double resolution=-1) |
kernel::Restraint * | create_excluded_volume_restraint (const Selections &s) |
Create an excluded volume restraint for a list of selections. More... | |
Hierarchy | create_fragment (const Hierarchies &ps) |
Create a fragment containing the specified nodes. More... | |
kernel::Restraint * | create_internal_connectivity_restraint (const Selection &s, double k, std::string name="Connectivity%1%") |
Create a restraint connecting the selection. More... | |
kernel::Restraint * | create_internal_connectivity_restraint (const Selection &s, double x0, double k, std::string name="Connectivity%1%") |
Create a restraint connecting the selection. More... | |
Hierarchy | create_protein (kernel::Model *m, std::string name, double target_radius, int number_of_residues, int first_residue_index=0, double volume=-1) |
Create a coarse grained molecule. More... | |
Hierarchy | create_protein (kernel::Model *m, std::string name, double target_radius, const Ints domain_boundaries) |
IMP::core::RigidBody | create_rigid_body (const Hierarchies &h, std::string name=std::string("created rigid body")) |
Rigidify a molecule or collection of molecules. More... | |
IMP::core::RigidBody | create_rigid_body (Hierarchy h) |
Hierarchy | create_simplified_assembly_from_volume (Hierarchy h, double resolution) |
Hierarchy | create_simplified_from_volume (Hierarchy h, double resolution) |
void | destroy (Hierarchy d) |
Delete the Hierarchy. More... | |
void | destroy_bond (Bond b) |
Destroy the bond connecting two particles. More... | |
CHARMMParameters * | get_all_atom_CHARMM_parameters () |
Atom | get_atom (Residue rd, AtomType at) |
Return a particle atom from the residue. More... | |
bool | get_atom_type_exists (std::string name) |
Return true if that atom type already exists. More... | |
Bond | get_bond (Bonded a, Bonded b) |
Get the bond between two particles. More... | |
algebra::BoundingBoxD< 3 > | get_bounding_box (const Hierarchy &h) |
Get a bounding box for the Hierarchy. More... | |
algebra::Sphere3D | get_bounding_sphere (const Hierarchy &h) |
Hierarchies | get_by_type (Hierarchy mhd, GetByType t) |
Chain | get_chain (Hierarchy h) |
std::string | get_chain_id (Hierarchy h) |
Atoms | get_charmm_untyped_atoms (Hierarchy hierarchy) |
Get all atoms in the Hierarchy that do not have CHARMM types. More... | |
FloatPair | get_component_placement_score (const core::XYZs &ref1, const core::XYZs &ref2, const core::XYZs &mdl1, const core::XYZs &mdl2) |
Measure the difference between two placements of the same set of points. More... | |
int | get_copy_index (Hierarchy h) |
Walk up the hierarchy to find the current copy index. More... | |
double | get_diffusion_angle (double D, double dtfs) |
double | get_diffusion_coefficient (const algebra::Vector3Ds &displacements, double dt) |
double | get_diffusion_coefficient_from_cm2_per_second (double din) |
double | get_diffusion_length (double D, double t) |
double | get_diffusion_length (double D, double force, double t, double temp=273) |
template<class Vector3DsOrXYZs0 , class Vector3DsOrXYZs1 > | |
double | get_drms (const Vector3DsOrXYZs0 &m1, const Vector3DsOrXYZs1 &m2) |
template<class Vector3DsOrXYZs0 , class Vector3DsOrXYZs1 > | |
double | get_drmsd (const Vector3DsOrXYZs0 &m0, const Vector3DsOrXYZs1 &m1) |
Calculate distance the root mean square deviation between two sets of 3D points. More... | |
template<class Vector3DsOrXYZs0 , class Vector3DsOrXYZs1 > | |
double | get_drmsd_Q (const Vector3DsOrXYZs0 &m0, const Vector3DsOrXYZs1 &m1, double threshold) |
double | get_einstein_diffusion_coefficient (double r) |
double | get_einstein_rotational_diffusion_coefficient (double r) |
Element | get_element_for_atom_type (AtomType at) |
ElementTable & | get_element_table () |
CHARMMParameters * | get_heavy_atom_CHARMM_parameters () |
HierarchyTree | get_hierarchy_tree (Hierarchy h) |
Get a graph for the passed Hierarchy. More... | |
Bonds | get_internal_bonds (Hierarchy mhd) |
Get the bonds internal to this tree. More... | |
bool | get_is_heterogen (Hierarchy h) |
Return true if the piece of hierarchy should be classified as a heterogen. More... | |
double | get_kd (double na, double nb, double nab, double volume) |
double | get_kt (double T) |
Return kT for a given temperature. More... | |
Hierarchies | get_leaves (const Selection &h) |
Hierarchies | get_leaves (Hierarchy h) |
Hierarchies | get_leaves (const Hierarchies &h) |
double | get_mass (ResidueType c) |
Get the mass from the residue type. More... | |
double | get_mass (const Selection &s) |
Get the total mass of a hierarchy, in Daltons. More... | |
double | get_maximum_time_step_estimate (BrownianDynamics *bd) |
double | get_molarity (double n, double volume) |
std::string | get_molecule_name (Hierarchy h) |
template<class Vector3DsOrXYZs0 , class Vector3DsOrXYZs1 > | |
double | get_native_overlap (const Vector3DsOrXYZs0 &m1, const Vector3DsOrXYZs1 &m2, double threshold) |
Computes the native overlap between two sets of 3D points. More... | |
Hierarchy | get_next_residue (Residue rd) |
Return the residue from the same chain with one higher index. More... | |
char | get_one_letter_code (ResidueType c) |
Get the 1-letter amino acid code from the residue type. More... | |
double | get_pairwise_rmsd_score (const core::XYZs &ref1, const core::XYZs &ref2, const core::XYZs &mdl1, const core::XYZs &mdl2) |
Measure the RMSD between two placements of the same set of points. More... | |
Atoms | get_phi_dihedral_atoms (Residue rd) |
Return the atoms comprising the phi dihedral. More... | |
FloatPair | get_placement_score (const core::XYZs &source, const core::XYZs &target) |
Measure the difference between two placements of the same set of points. More... | |
Hierarchy | get_previous_residue (Residue rd) |
Return the residue from the same chain with one lower index. More... | |
Atoms | get_psi_dihedral_atoms (Residue rd) |
Return the atoms comprising the psi dihedral. More... | |
double | get_radius_of_gyration (const Selection &s) |
double | get_radius_of_gyration (const kernel::ParticlesTemp &ps) |
Residue | get_residue (Atom d, bool nothrow=false) |
Return the Residue containing this atom. More... | |
Hierarchy | get_residue (Hierarchy mhd, unsigned int index) |
Get the residue with the specified index. More... | |
ResidueType | get_residue_type (char c) |
Get the residue type from the 1-letter amino acid code. More... | |
double | get_resolution (kernel::Model *m, kernel::ParticleIndex pi) |
double | get_resolution (Hierarchy h) |
template<class Vector3DsOrXYZs0 , class Vector3DsOrXYZs1 > | |
double | get_rigid_bodies_drms (const Vector3DsOrXYZs0 &m1, const Vector3DsOrXYZs1 &m2, const IMP::IntRanges &ranges) |
double | get_rmsd (const core::XYZs &s0, const core::XYZs &s1) |
template<class Vector3DsOrXYZs0 , class Vector3DsOrXYZs1 > | |
double | get_rmsd (const Vector3DsOrXYZs0 &m1, const Vector3DsOrXYZs1 &m2, const algebra::Transformation3D &tr_for_second=algebra::get_identity_transformation_3d()) |
double | get_rmsd (const Selection &s0, const Selection &s1) |
double | get_rmsd (const Selection &s0, const Selection &s1, const algebra::Transformation3D &tr_for_second) |
double | get_rmsd (const algebra::Vector3Ds &s0, const algebra::Vector3Ds &s1, const algebra::Transformation3D &tr_for_second) |
double | get_rmsd (const core::XYZs &s0, const core::XYZs &s1, const IMP::algebra::Transformation3D &tr_for_second) |
double | get_rmsd_transforming_first (const IMP::algebra::Transformation3D &tr, const core::XYZs &s0, const core::XYZs &s1) |
double | get_rmsd_transforming_first (const algebra::Transformation3D &tr, const Selection &s0, const Selection &s1) |
Hierarchy | get_root (Hierarchy h) |
Return the root of the hierarchy. More... | |
double | get_rotational_diffusion_coefficient (const algebra::Rotation3Ds &displacements, double dt) |
Float | get_secondary_structure_match_score (SecondaryStructureResidue ssr1, SecondaryStructureResidue ssr2) |
Compares the secondary structure probabilities of two residues. More... | |
int | get_state_index (Hierarchy h) |
Walk up the hierarchy to find the current state. More... | |
double | get_surface_area (const Selection &s) |
Get the total surface area of a hierarchy, in square angstroms. More... | |
HierarchyTreeVertexIndex | get_vertex_index (const HierarchyTree &g) |
double | get_volume (const Selection &s) |
Get the total volume of a hierarchy, in cubic angstroms. More... | |
void | remove_charmm_untyped_atoms (Hierarchy hierarchy) |
Remove any atom from the Hierarchy that does not have a CHARMM type. More... | |
void | setup_as_approximation (kernel::Particle *h, const kernel::ParticlesTemp &other) |
void | setup_as_approximation (Hierarchy h) |
SecondaryStructureResidue | setup_coarse_secondary_structure_residue (const kernel::Particles &ssr_ps, kernel::Model *mdl, bool winner_takes_all_per_res=false) |
Coarsen some SecondaryStructureResidues. More... | |
SecondaryStructureResidues | setup_coarse_secondary_structure_residues (const kernel::Particles &ssr_ps, kernel::Model *mdl, int coarse_factor, int start_res_num, bool winner_takes_all_per_res=false) |
void | show (Hierarchy h, std::ostream &out=std::cout) |
Print out a molecular hierarchy. More... | |
void | show_as_graphviz (const HierarchyTree &name, base::TextOutput out) |
void | transform (Hierarchy h, const algebra::Transformation3D &tr) |
Transform a hierarchy. This is aware of rigid bodies. More... | |
Variables | |
const double | ALL_RESOLUTIONS |
const AtomType | AT_CA |
const AtomType | AT_N |
const AtomType | AT_UNKNOWN |
const ResidueType | GLY |
const ResidueType | UNK |
Unknown residue. More... | |
Standard module functions | |
All | |
std::string | get_module_version () |
std::string | get_module_name () |
std::string | get_data_path (std::string file_name) |
Return the full path to installed data. More... | |
std::string | get_example_path (std::string file_name) |
Return the path to installed example data for this module. More... | |
Dope scoring | |
| |
typedef score_functor::DopeType | DopeType |
void | add_dope_score_data (atom::Hierarchy h) |
Estimator Functions | |
These functions allow you to estimate physical quantities relating to biomolecules. | |
enum | ProteinDensityReference { ALBER, HARPAZ, ANDERSSON, TSAI, QUILLIN, SQUIRE } |
double | get_protein_density_from_reference (ProteinDensityReference densityReference) |
double | get_volume_from_mass (double m, ProteinDensityReference ref=ALBER) |
Estimate the volume of a protein from its mass. More... | |
double | get_mass_from_volume (double v, ProteinDensityReference ref=ALBER) |
Estimate the mass of a protein from its volume. More... | |
double | get_mass_from_number_of_residues (unsigned int num_aa) |
Estimate the mass of a protein from the number of amino acids. More... | |
double | get_volume_from_residue_type (ResidueType rt) |
Return an estimate for the volume of a given residue. More... | |
Simplification along backbone | |
These two methods create a simplified version of a molecule by merging residues sequentially. In one case every n residues are merged, in the other, the intervals are passed manually. The resulting molecule is not optimized by default and has no restraints automatically created. At the moment, the calls only support unmodified hierarchies loaded by read_pdb() which have only protein or DNA members. They return Hierarchy() if the input chain is empty. If keep_detailed is true, then the original high resolution structure particles are added as children of the simplified structure. | |
Hierarchy | create_simplified_along_backbone (Hierarchy input, int num_res, bool keep_detailed=false) |
Hierarchy | create_simplified_along_backbone (Chain input, const IntRanges &residue_segments, bool keep_detailed=false) |
Finding information | |
Get the attribute of the given particle or throw a ValueException if it is not applicable. The particle with the given information must be above the passed node. | |
Ints | get_residue_indexes (Hierarchy h) |
ResidueType | get_residue_type (Hierarchy h) |
char | get_chain_id_char (Hierarchy h) |
AtomType | get_atom_type (Hierarchy h) |
std::string | get_domain_name (Hierarchy h) |
Mol2 IO | |
The read function produces a hierarchy containing the molecule. The write hierarchy writes all the Residue types in the hierarchy to the file. | |
Hierarchy | read_mol2 (base::TextInput mol2_file, kernel::Model *model, Mol2Selector *mol2sel=nullptr) |
Create a hierarchy from a Mol2 file. More... | |
void | write_mol2 (Hierarchy rhd, base::TextOutput file_name) |
Write a ligand hierarchy as a mol2 file. More... | |
PDB Reading | |
The read PDB methods produce a hierarchy that looks as follows:
Waters are currently dropped if they are ATOM records. This can be fixed. The read_pdb() functions should successfully parse all valid PDB files. It can produce warnings on files which are not valid. It will attempt to read such files, but all bets are off. When reading PDBs, PDBSelector objects can be used to choose to only process certain record types. See the class documentation for more information. When no PDB selector is supplied for reading, the NonWaterPDBSelector is used. Set the IMP::LogLevel to VERBOSE to see details of parse errors. | |
PDBSelector * | get_default_pdb_selector () |
Hierarchy | read_pdb (base::TextInput input, kernel::Model *model, PDBSelector *selector=get_default_pdb_selector(), bool select_first_model=true) |
void | read_pdb (base::TextInput input, int model, Hierarchy h) |
Hierarchies | read_multimodel_pdb (base::TextInput input, kernel::Model *model, PDBSelector *selector=get_default_pdb_selector()) |
PDB Writing | |
The methods to write a PDB expects a Hierarchy that looks as follows: All Residue particles that have a Chain particle as an ancestor are considered part of a protein, DNA or RNA, ones without are considered heterogens. The functions produce files that are not valid PDB files, eg only ATOM/HETATM lines are printed for all Atom particles in the hierarchy. Complain if your favorite program can't read them and we might fix it. | |
void | write_pdb (const Selection &mhd, base::TextOutput out, unsigned int model=1) |
void | write_pdb_of_c_alphas (const Selection &mhd, base::TextOutput out, unsigned int model=1) |
Write a hierarchy to a PDB as C_alpha atoms. More... | |
void | write_multimodel_pdb (const Hierarchies &mhd, base::TextOutput out) |
Protein-ligand scoring | |
| |
typedef Key< 783462, false > | ProteinLigandType |
typedef IMP::base::Vector < ProteinLigandType > | ProteinLigandTypes |
typedef IMP::base::Vector < IMP::base::Pointer < ProteinLigandAtomPairScore > > | ProteinLigandAtomPairScores |
typedef IMP::base::Vector < IMP::base::WeakPointer < ProteinLigandAtomPairScore > > | ProteinLigandAtomPairScoresTemp |
typedef IMP::base::Vector < IMP::base::Pointer < ProteinLigandRestraint > > | ProteinLigandRestraints |
typedef IMP::base::Vector < IMP::base::WeakPointer < ProteinLigandRestraint > > | ProteinLigandRestraintsTemp |
void | add_protein_ligand_score_data (Hierarchy h) |
Python Only | |
The following functions are only available in Python as the equivalent C++ functionality is provided via template functions or in other ways that don't directly map to Python. | |
void | show_molecular_hierarchy (Hierarchy h) |
Print out the molecular hierarchy. More... | |
PSIPRED reading | |
Reads in PSIPRED results and decorates particles as SecondaryStructureResidues. Currently assuming order of ps matches file. | |
SecondaryStructureResidues | read_psipred (base::TextInput inf, kernel::Model *mdl) |
SecondaryStructureResidues | read_psipred (base::TextInput inf, kernel::Particles ps) |
typedef IMP::base::Vector<IMP::base::Pointer< AngleSingletonScore > > IMP::atom::AngleSingletonScores |
Store a set of objects.
Definition at line 45 of file AngleSingletonScore.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< AngleSingletonScore > > IMP::atom::AngleSingletonScoresTemp |
Pass a set of objects.
Definition at line 45 of file AngleSingletonScore.h.
typedef IMP::base::Vector< AtomType > IMP::atom::AtomTypes |
Pass or store a set of AtomType .
Definition at line 29 of file atom/Atom.h.
typedef IMP::base::Vector<IMP::base::Pointer< BerendsenThermostatOptimizerState > > IMP::atom::BerendsenThermostatOptimizerStates |
Store a set of objects.
Definition at line 70 of file BerendsenThermostatOptimizerState.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< BerendsenThermostatOptimizerState > > IMP::atom::BerendsenThermostatOptimizerStatesTemp |
Pass a set of objects.
Definition at line 70 of file BerendsenThermostatOptimizerState.h.
Store a set of objects.
Definition at line 38 of file BondedPairFilter.h.
typedef IMP::base::Vector<IMP::base::Pointer< BondEndpointsRefiner > > IMP::atom::BondEndpointsRefiners |
Store a set of objects.
Definition at line 37 of file BondEndpointsRefiner.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< BondEndpointsRefiner > > IMP::atom::BondEndpointsRefinersTemp |
Pass a set of objects.
Definition at line 37 of file BondEndpointsRefiner.h.
Store a set of objects.
Definition at line 61 of file BondPairContainer.h.
Store a set of objects.
Definition at line 42 of file BondSingletonScore.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< BondSingletonScore > > IMP::atom::BondSingletonScoresTemp |
Pass a set of objects.
Definition at line 42 of file BondSingletonScore.h.
typedef IMP::base::Vector< CHARMMAngle > IMP::atom::CHARMMAngles |
Pass or store a set of CHARMMAngle .
Definition at line 169 of file charmm_topology.h.
Pass or store a set of CHARMMAtomTopology .
Definition at line 50 of file charmm_topology.h.
Pass or store a set of CHARMMBondEndpoint .
Definition at line 90 of file charmm_topology.h.
Pass or store a set of CHARMMBondParameters .
Definition at line 38 of file CHARMMParameters.h.
typedef IMP::base::Vector< CHARMMBond > IMP::atom::CHARMMBonds |
Pass or store a set of CHARMMBond .
Definition at line 167 of file charmm_topology.h.
Pass or store a set of CHARMMDihedralParameters .
Definition at line 52 of file CHARMMParameters.h.
Pass or store a set of CHARMMDihedral .
Definition at line 171 of file charmm_topology.h.
typedef IMP::base::Vector<IMP::base::Pointer< CHARMMIdealResidueTopology > > IMP::atom::CHARMMIdealResidueTopologies |
Store a set of objects.
Definition at line 358 of file charmm_topology.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< CHARMMIdealResidueTopology > > IMP::atom::CHARMMIdealResidueTopologiesTemp |
Pass a set of objects.
Definition at line 358 of file charmm_topology.h.
Pass or store a set of CHARMMInternalCoordinate .
Definition at line 239 of file charmm_topology.h.
Store a set of objects.
Definition at line 345 of file CHARMMParameters.h.
Store a set of objects.
Definition at line 413 of file charmm_topology.h.
typedef IMP::base::Vector<IMP::base::Pointer< CHARMMResidueTopology > > IMP::atom::CHARMMResidueTopologies |
Store a set of objects.
Definition at line 465 of file charmm_topology.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< CHARMMResidueTopology > > IMP::atom::CHARMMResidueTopologiesTemp |
Pass a set of objects.
Definition at line 465 of file charmm_topology.h.
typedef IMP::base::Vector<IMP::base::Pointer< CHARMMResidueTopologyBase > > IMP::atom::CHARMMResidueTopologyBases |
Store a set of objects.
Definition at line 324 of file charmm_topology.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< CHARMMResidueTopologyBase > > IMP::atom::CHARMMResidueTopologyBasesTemp |
Pass a set of objects.
Definition at line 324 of file charmm_topology.h.
typedef IMP::base::Vector<IMP::base::Pointer< CHARMMSegmentTopology > > IMP::atom::CHARMMSegmentTopologies |
Store a set of objects.
Definition at line 48 of file charmm_segment_topology.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< CHARMMSegmentTopology > > IMP::atom::CHARMMSegmentTopologiesTemp |
Pass a set of objects.
Definition at line 48 of file charmm_segment_topology.h.
Store a set of objects.
Definition at line 244 of file charmm_segment_topology.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< CHARMMTopology > > IMP::atom::CHARMMTopologiesTemp |
Pass a set of objects.
Definition at line 244 of file charmm_segment_topology.h.
Store a set of objects.
Definition at line 58 of file CoulombPairScore.h.
Store a set of objects.
Definition at line 39 of file CoverBond.h.
typedef IMP::base::Vector<IMP::base::Pointer< DihedralSingletonScore > > IMP::atom::DihedralSingletonScores |
Store a set of objects.
Definition at line 42 of file DihedralSingletonScore.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< DihedralSingletonScore > > IMP::atom::DihedralSingletonScoresTemp |
Pass a set of objects.
Definition at line 42 of file DihedralSingletonScore.h.
typedef boost::graph IMP::atom::HierarchyTree |
A graph for representing a Hierarchy so you can view it nicely.
See Graphs in IMP for more information.
Definition at line 170 of file hierarchy_tools.h.
typedef IMP::base::Vector<IMP::base::Pointer< ImproperSingletonScore > > IMP::atom::ImproperSingletonScores |
Store a set of objects.
Definition at line 50 of file ImproperSingletonScore.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< ImproperSingletonScore > > IMP::atom::ImproperSingletonScoresTemp |
Pass a set of objects.
Definition at line 50 of file ImproperSingletonScore.h.
typedef IMP::base::Vector<IMP::base::Pointer< LangevinThermostatOptimizerState > > IMP::atom::LangevinThermostatOptimizerStates |
Store a set of objects.
Definition at line 59 of file LangevinThermostatOptimizerState.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< LangevinThermostatOptimizerState > > IMP::atom::LangevinThermostatOptimizerStatesTemp |
Pass a set of objects.
Definition at line 59 of file LangevinThermostatOptimizerState.h.
typedef IMP::base::Vector<IMP::base::Pointer< LennardJonesPairScore > > IMP::atom::LennardJonesPairScores |
Store a set of objects.
Definition at line 89 of file LennardJonesPairScore.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< LennardJonesPairScore > > IMP::atom::LennardJonesPairScoresTemp |
Pass a set of objects.
Definition at line 89 of file LennardJonesPairScore.h.
typedef IMP::base::Vector<IMP::base::Pointer< ProteinLigandAtomPairScore > > IMP::atom::ProteinLigandAtomPairScores |
Store a set of objects.
Definition at line 65 of file protein_ligand_score.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< ProteinLigandAtomPairScore > > IMP::atom::ProteinLigandAtomPairScoresTemp |
Pass a set of objects.
Definition at line 65 of file protein_ligand_score.h.
typedef IMP::base::Vector<IMP::base::Pointer< ProteinLigandRestraint > > IMP::atom::ProteinLigandRestraints |
Store a set of objects.
Definition at line 80 of file protein_ligand_score.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< ProteinLigandRestraint > > IMP::atom::ProteinLigandRestraintsTemp |
Pass a set of objects.
Definition at line 80 of file protein_ligand_score.h.
typedef Key<783462, false> IMP::atom::ProteinLigandType |
The marker to identify the atom types.
Definition at line 34 of file protein_ligand_score.h.
Pass or store a set of ProteinLigandType .
Definition at line 39 of file protein_ligand_score.h.
typedef IMP::base::Vector<IMP::base::Pointer< RemoveRigidMotionOptimizerState > > IMP::atom::RemoveRigidMotionOptimizerStates |
Store a set of objects.
Definition at line 45 of file RemoveRigidMotionOptimizerState.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< RemoveRigidMotionOptimizerState > > IMP::atom::RemoveRigidMotionOptimizerStatesTemp |
Pass a set of objects.
Definition at line 45 of file RemoveRigidMotionOptimizerState.h.
typedef IMP::base::Vector<IMP::base::Pointer< RemoveTranslationOptimizerState > > IMP::atom::RemoveTranslationOptimizerStates |
Store a set of objects.
Definition at line 43 of file RemoveTranslationOptimizerState.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< RemoveTranslationOptimizerState > > IMP::atom::RemoveTranslationOptimizerStatesTemp |
Pass a set of objects.
Definition at line 43 of file RemoveTranslationOptimizerState.h.
typedef IMP::base::Vector< ResidueType > IMP::atom::ResidueTypes |
Pass or store a set of ResidueType .
typedef IMP::base::Vector< Selection > IMP::atom::Selections |
Pass or store a set of Selection .
Definition at line 217 of file Selection.h.
Store a set of objects.
Definition at line 181 of file Simulator.h.
typedef IMP::base::Vector<IMP::base::Pointer< VelocityScalingOptimizerState > > IMP::atom::VelocityScalingOptimizerStates |
Store a set of objects.
Definition at line 55 of file VelocityScalingOptimizerState.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< VelocityScalingOptimizerState > > IMP::atom::VelocityScalingOptimizerStatesTemp |
Pass a set of objects.
Definition at line 55 of file VelocityScalingOptimizerState.h.
enum IMP::atom::Element |
enum IMP::atom::GetByType |
The different types which can be passed to get_by_type()
Definition at line 367 of file atom/Hierarchy.h.
Several protein density value references that have been proposed in the literature.
Definition at line 33 of file estimates.h.
Eventually, other types of representation will be supported, eg Gaussians or density maps.
Definition at line 26 of file Representation.h.
AtomType IMP::atom::add_atom_type | ( | std::string | name, |
Element | e | ||
) |
void IMP::atom::add_bonds | ( | Hierarchy | d, |
const ForceFieldParameters * | ffp = get_all_atom_CHARMM_parameters() |
||
) |
Add bonds using definitions from given force field parameters. Note that, at the moment, all added bonds are reported as IMP::Bond::SINGLE, whether or not they actually are.
void IMP::atom::add_dope_score_data | ( | atom::Hierarchy | h | ) |
Add the dope atom types to the atoms in the hierarchy.
void IMP::atom::add_protein_ligand_score_data | ( | Hierarchy | h | ) |
Add the data needed to use ProteinLigandAtomPairScore with the passed Hierarchy.
void IMP::atom::add_radii | ( | Hierarchy | d, |
const ForceFieldParameters * | ffp = get_all_atom_CHARMM_parameters() , |
||
FloatKey | radius_key = FloatKey("radius") |
||
) |
Add vdW radius from given force field.
Hierarchy IMP::atom::create_clone | ( | Hierarchy | d | ) |
Hierarchy IMP::atom::create_clone_one | ( | Hierarchy | d | ) |
IMP::core::RigidBody IMP::atom::create_compatible_rigid_body | ( | Hierarchy | h, |
Hierarchy | reference | ||
) |
Rigidify a molecule or collection of molecules.
This method is identical to create_rigid_body() except that the chosen reference frame is aligned with that of reference (which must have exactly the same set of particles). This allows one to make sure the rigid body is equivalent when you have several copies of the same molecule.
kernel::Restraint* IMP::atom::create_connectivity_restraint | ( | const Selections & | s, |
double | k, | ||
std::string | name = "Connectivity%1%" |
||
) |
Create a restraint connecting the selections.
If one or more of the selections is a rigid body, this will be used to accelerate the computation.
kernel::Restraint* IMP::atom::create_connectivity_restraint | ( | const Selections & | s, |
double | x0, | ||
double | k, | ||
std::string | name = "Connectivity%1%" |
||
) |
Create a restraint connecting the selections.
The particles are allowed to be apart by x0 and still count as connected.
If one or more of the selections is a rigid body, this will be used to accelerate the computation.
core::XYZR IMP::atom::create_cover | ( | const Selection & | s, |
std::string | name = std::string() |
||
) |
Create an XYZR particle which always includes the particles in the selection in its bounding volume. If all the particles in the selection are part of the same rigid body, then the created particle is added as part of that rigid body. Otherwise it uses an IMP::core::Cover to maintain the cover property.
Doing this can be a useful way to accelerate computations when it is OK to replace a potentially complicated set of geometries represented by the selection with a much simpler one.
kernel::Restraint* IMP::atom::create_distance_restraint | ( | const Selection & | n0, |
const Selection & | n1, | ||
double | x0, | ||
double | k, | ||
std::string | name = "Distance%1%" |
||
) |
Create a distance restraint between the selections.
This restraint applies a harmonic to the minimum distance between a particle in selection n0 and a particle in selection n1.
If one or more of the selections is a rigid body, this will be used to accelerate the computation.
kernel::Restraint* IMP::atom::create_excluded_volume_restraint | ( | const Hierarchies & | hs, |
double | resolution = -1 |
||
) |
Create an excluded volume restraint for the included molecules. If a value is provided for resolution, then something less than the full resolution representation will be used.
If one or more of the selections is a rigid body, this will be used to accelerate the computation.
kernel::Restraint* IMP::atom::create_excluded_volume_restraint | ( | const Selections & | s | ) |
Create an excluded volume restraint for a list of selections.
Hierarchy IMP::atom::create_fragment | ( | const Hierarchies & | ps | ) |
Create a fragment containing the specified nodes.
A particle representing the fragment is created and initialized.
The Fragment is inserted as a child of the parent (and the particles are removed). The particles become children of the fragment.
ValueException | If all the particles do not have the same parent. |
kernel::Restraint* IMP::atom::create_internal_connectivity_restraint | ( | const Selection & | s, |
double | k, | ||
std::string | name = "Connectivity%1%" |
||
) |
Create a restraint connecting the selection.
If one or more of the selections is a rigid body, this will be used to accelerate the computation.
kernel::Restraint* IMP::atom::create_internal_connectivity_restraint | ( | const Selection & | s, |
double | x0, | ||
double | k, | ||
std::string | name = "Connectivity%1%" |
||
) |
Create a restraint connecting the selection.
The particles are allowed to be apart by x0 and still count as connected.
If one or more of the selections is a rigid body, this will be used to accelerate the computation.
Hierarchy IMP::atom::create_protein | ( | kernel::Model * | m, |
std::string | name, | ||
double | target_radius, | ||
int | number_of_residues, | ||
int | first_residue_index = 0 , |
||
double | volume = -1 |
||
) |
Create a coarse grained molecule.
The coarse grained model is created with a number of spheres based on the resolution and the volume. If the volume is not provided it is estimated based on the number of residues. The protein is created as a molecular hierarchy rooted at p. The leaves are Domain particles with appropriate residue indexes stored and are XYZR particles.
Volume is, as usual, in cubic angstroms.
Currently the function creates a set of balls with radii no greater than target_radius which overlap by 20% and have a volume of their union equal to the passed volume.
The coordinates of the balls defining the protein are not optimized by default, and have garbage coordinate values.
Hierarchy IMP::atom::create_protein | ( | kernel::Model * | m, |
std::string | name, | ||
double | target_radius, | ||
const Ints | domain_boundaries | ||
) |
Like the former create_protein(), but it enforces domain splits at the provided domain boundaries. The domain boundaries should be the start of the first domain, any boundaries, and then one past the end of the last domain.
IMP::core::RigidBody IMP::atom::create_rigid_body | ( | const Hierarchies & | h, |
std::string | name = std::string("created rigid body") |
||
) |
Rigidify a molecule or collection of molecules.
The rigid body created has all the leaves as members and a member rigid body for each internal node in the tree. The particle created to be the rigid body is returned.
A name can be passed as it is not easy to automatically pick a decent name.
IMP::core::RigidBody IMP::atom::create_rigid_body | ( | Hierarchy | h | ) |
Hierarchy IMP::atom::create_simplified_along_backbone | ( | Hierarchy | input, |
int | num_res, | ||
bool | keep_detailed = false |
||
) |
Simplify every num_res into one particle.
Hierarchy IMP::atom::create_simplified_along_backbone | ( | Chain | input, |
const IntRanges & | residue_segments, | ||
bool | keep_detailed = false |
||
) |
Simplify by breaking at the boundaries provided.
Hierarchy IMP::atom::create_simplified_assembly_from_volume | ( | Hierarchy | h, |
double | resolution | ||
) |
Create a new hierarchy that approximates the volume occupied by the old one.
This function is like create_simplified_from_volume() except that the result is divided into Chain and Molecule bits. It assumes that all geometry is rooted under a chain or a molecule.
Hierarchy IMP::atom::create_simplified_from_volume | ( | Hierarchy | h, |
double | resolution | ||
) |
Create a new hierarchy that approximates the volume occupied by the old one.
The new hierarchy will contain a number of Fragment particles whose surface (and hence volume) approximates the input hierarchy in the sense of IMP::algebra::get_simplified_from_volume().
The resulting representation has approximately the desired resolution in the sense of get_resolution().
Residue indexes are assigned to the created Fragments from the original residue indexes, and particles containing adjacent residues are connected by bonds.
void IMP::atom::destroy | ( | Hierarchy | d | ) |
Delete the Hierarchy.
All bonds connecting to these atoms are destroyed as are hierarchy links in the Hierarchy and the particles are removed from the kernel::Model. If this particle has a parent, it is removed from the parent.
CHARMMParameters* IMP::atom::get_all_atom_CHARMM_parameters | ( | ) |
The default CHARMM parameters support normal amino acid and nucleic acid residues and the atoms found in them. To use CHARMM with heterogens or non-standard residues, a different CHARMM parameters file must be used.
Atom IMP::atom::get_atom | ( | Residue | rd, |
AtomType | at | ||
) |
bool IMP::atom::get_atom_type_exists | ( | std::string | name | ) |
Return true if that atom type already exists.
algebra::BoundingBoxD<3> IMP::atom::get_bounding_box | ( | const Hierarchy & | h | ) |
Get a bounding box for the Hierarchy.
This bounding box is that of the highest (in the CS sense of a tree growing down from the root) cut through the tree where each node in the cut has x,y,z, and r. That is, if the root has x,y,z,r then it is the bounding box of that sphere. If only the leaves have radii, it is the bounding box of the leaves. If no such cut exists, the behavior is undefined.
algebra::Sphere3D IMP::atom::get_bounding_sphere | ( | const Hierarchy & | h | ) |
See get_bounding_box() for more details.
Hierarchies IMP::atom::get_by_type | ( | Hierarchy | mhd, |
GetByType | t | ||
) |
Gather all the molecular particles of a certain level in the molecular hierarchy.
Chain IMP::atom::get_chain | ( | Hierarchy | h | ) |
Get the containing chain or Chain() if there is none
std::string IMP::atom::get_chain_id | ( | Hierarchy | h | ) |
Walk up the hierarchy to determine the chain id.
char IMP::atom::get_chain_id_char | ( | Hierarchy | h | ) |
Definition at line 124 of file hierarchy_tools.h.
Atoms IMP::atom::get_charmm_untyped_atoms | ( | Hierarchy | hierarchy | ) |
Get all atoms in the Hierarchy that do not have CHARMM types.
FloatPair IMP::atom::get_component_placement_score | ( | const core::XYZs & | ref1, |
const core::XYZs & | ref2, | ||
const core::XYZs & | mdl1, | ||
const core::XYZs & | mdl2 | ||
) |
Measure the difference between two placements of the same set of points.
[in] | ref1 | The reference placement of the first component represented by XYZ coordinates |
[in] | ref2 | The reference placement of the second component represented by XYZ coordinates |
[in] | mdl1 | The modeled placement of the first component represented by XYZ coordinates |
[in] | mdl2 | The modeled placement of the second component represented by XYZ coordinates |
int IMP::atom::get_copy_index | ( | Hierarchy | h | ) |
Walk up the hierarchy to find the current copy index.
std::string IMP::atom::get_data_path | ( | std::string | file_name | ) |
Return the full path to installed data.
Each module has its own data directory, so be sure to use the version of this function in the correct module. To read the data file "data_library" that was placed in the data
directory of module "mymodule", do something like
This will ensure that the code works when IMP
is installed or used via the setup_environment.sh
script.
double IMP::atom::get_diffusion_angle | ( | double | D, |
double | dtfs | ||
) |
Get the standard deviation of the diffusion angle change given the rigid body diffusion coefficient and the time step dtfs.
double IMP::atom::get_diffusion_coefficient | ( | const algebra::Vector3Ds & | displacements, |
double | dt | ||
) |
Estimate the diffusion coefficient of a particle from a list of displacements each taken after the given time step dt.
double IMP::atom::get_diffusion_length | ( | double | D, |
double | t | ||
) |
Return the standard deviation for the Brownian step given the diffusion coefficient D and the time step t.
double IMP::atom::get_diffusion_length | ( | double | D, |
double | force, | ||
double | t, | ||
double | temp = 273 |
||
) |
Return the scale for diffusion under the specified force, the diffusion coefficient D and the time step t.
D | diffusion coefficient in Angstrom^2/femtosecond |
force | applied force |
t | time step in femtoseconds |
temp | temperature in Kelvin |
double IMP::atom::get_drms | ( | const Vector3DsOrXYZs0 & | m1, |
const Vector3DsOrXYZs1 & | m2 | ||
) |
Distance-RMS between two sets of points, defined as: sqrt( sum[ (d1ij**2 - d2ij**2)**2]/(4 * sum[d1ij**2]) ) (Levitt, 1992) d1ij - distance between points i and j in set m1 (the "reference" set) d2ij - distance between points i and j in set m2
[in] | m1 | set of points |
[in] | m2 | set of points |
Definition at line 186 of file atom/distance.h.
double IMP::atom::get_drmsd | ( | const Vector3DsOrXYZs0 & | m0, |
const Vector3DsOrXYZs1 & | m1 | ||
) |
Calculate distance the root mean square deviation between two sets of 3D points.
See generic geometry for more information.
Definition at line 88 of file atom/distance.h.
double IMP::atom::get_drmsd_Q | ( | const Vector3DsOrXYZs0 & | m0, |
const Vector3DsOrXYZs1 & | m1, | ||
double | threshold | ||
) |
Calculate the distance root mean square deviation between two sets of 3D points with a distance threshold
See generic geometry for more information.
Definition at line 124 of file atom/distance.h.
double IMP::atom::get_einstein_diffusion_coefficient | ( | double | r | ) |
Return the prediction diffusion coefficient in Angstrom squared per femtosecond given a radius in angstrom. See wikipedia for a reference and Wikipedia on Viscosity for the values of the viscosity of water used.
double IMP::atom::get_einstein_rotational_diffusion_coefficient | ( | double | r | ) |
Return the prediction diffusion coefficient in radians squared per femtosecond given a radius in angstrom.
std::string IMP::atom::get_example_path | ( | std::string | file_name | ) |
Return the path to installed example data for this module.
Each module has its own example directory, so be sure to use the version of this function in the correct module. For example to read the file example_protein.pdb
located in the examples
directory of the IMP::atom module, do
This will ensure that the code works when IMP
is installed or used via the setup_environment.sh
script.
CHARMMParameters* IMP::atom::get_heavy_atom_CHARMM_parameters | ( | ) |
The default CHARMM parameters support normal amino acid and nucleic acid residues and the atoms found in them. To use CHARMM with heterogens or non-standard residues, a different CHARMM parameters file must be used.
No hydrogen parameters are read.
HierarchyTree IMP::atom::get_hierarchy_tree | ( | Hierarchy | h | ) |
Bonds IMP::atom::get_internal_bonds | ( | Hierarchy | mhd | ) |
bool IMP::atom::get_is_heterogen | ( | Hierarchy | h | ) |
Return true if the piece of hierarchy should be classified as a heterogen.
For the purposes of classification, a heterogen is anything that
double IMP::atom::get_kd | ( | double | na, |
double | nb, | ||
double | nab, | ||
double | volume | ||
) |
Compute the concentration in molarity from the passed values
Definition at line 94 of file estimates.h.
double IMP::atom::get_kt | ( | double | T | ) |
Return kT for a given temperature.
Value taken from Wikipedia.
Definition at line 18 of file atom/constants.h.
Hierarchies IMP::atom::get_leaves | ( | const Selection & | h | ) |
Hierarchies IMP::atom::get_leaves | ( | Hierarchy | h | ) |
Definition at line 433 of file atom/Hierarchy.h.
Hierarchies IMP::atom::get_leaves | ( | const Hierarchies & | h | ) |
Definition at line 438 of file atom/Hierarchy.h.
double IMP::atom::get_mass | ( | ResidueType | c | ) |
Get the mass from the residue type.
double IMP::atom::get_mass | ( | const Selection & | s | ) |
Get the total mass of a hierarchy, in Daltons.
double IMP::atom::get_mass_from_number_of_residues | ( | unsigned int | num_aa | ) |
Estimate the mass of a protein from the number of amino acids.
We use an estimate of 110 Daltons per residue, following Chimera.
The mass is in Daltons.
double IMP::atom::get_mass_from_volume | ( | double | v, |
ProteinDensityReference | ref = ALBER |
||
) |
Estimate the mass of a protein from its volume.
[in] | v | the volume for which we want to output the corresponding mass |
[in] | ref | the protein density reference used in the computation. As a default ref is the estimate published in Alber et. al, Structure 2005. |
double IMP::atom::get_maximum_time_step_estimate | ( | BrownianDynamics * | bd | ) |
Repeatedly run the current model with Brownian dynamics at different time steps to try to find the maximum time step that can be used without the model exploding.
double IMP::atom::get_molarity | ( | double | n, |
double | volume | ||
) |
Compute the concentration in molarity from the passed values
Definition at line 87 of file estimates.h.
std::string IMP::atom::get_molecule_name | ( | Hierarchy | h | ) |
Walk up the hierarchy to determine the molecule name.
double IMP::atom::get_native_overlap | ( | const Vector3DsOrXYZs0 & | m1, |
const Vector3DsOrXYZs1 & | m2, | ||
double | threshold | ||
) |
Computes the native overlap between two sets of 3D points.
[in] | m1 | first set |
[in] | m2 | second set |
[in] | threshold | threshold distance (angstroms) for the calculation |
Definition at line 163 of file atom/distance.h.
Hierarchy IMP::atom::get_next_residue | ( | Residue | rd | ) |
Return the residue from the same chain with one higher index.
If no such residue exists, return Hierarchy().
The return type is Hierarchy since the particle representing the next residue might not be a Residue particle.
char IMP::atom::get_one_letter_code | ( | ResidueType | c | ) |
Get the 1-letter amino acid code from the residue type.
double IMP::atom::get_pairwise_rmsd_score | ( | const core::XYZs & | ref1, |
const core::XYZs & | ref2, | ||
const core::XYZs & | mdl1, | ||
const core::XYZs & | mdl2 | ||
) |
Measure the RMSD between two placements of the same set of points.
[in] | ref1 | The reference placement of the first component represented by XYZ coordinates |
[in] | ref2 | The reference placement of the second component represented by XYZ coordinates |
[in] | mdl1 | The modeled placement of the first component represented by XYZ coordinates |
[in] | mdl2 | The modeled placement of the second component represented by XYZ coordinates |
Atoms IMP::atom::get_phi_dihedral_atoms | ( | Residue | rd | ) |
Return the atoms comprising the phi dihedral.
If all atoms cannot be found, an empty list is returned.
FloatPair IMP::atom::get_placement_score | ( | const core::XYZs & | source, |
const core::XYZs & | target | ||
) |
Measure the difference between two placements of the same set of points.
[in] | source | The reference placement represented by XYZ coordinates |
[in] | target | The modeled placement represented by XYZ coordinates |
Hierarchy IMP::atom::get_previous_residue | ( | Residue | rd | ) |
Return the residue from the same chain with one lower index.
If no such residue exists, return Hierarchy().
double IMP::atom::get_protein_density_from_reference | ( | ProteinDensityReference | densityReference | ) |
returns the protein density value (in Da/A^3) associated with a given reference
Atoms IMP::atom::get_psi_dihedral_atoms | ( | Residue | rd | ) |
Return the atoms comprising the psi dihedral.
If all atoms cannot be found, an empty list is returned.
double IMP::atom::get_radius_of_gyration | ( | const Selection & | s | ) |
double IMP::atom::get_radius_of_gyration | ( | const kernel::ParticlesTemp & | ps | ) |
Compute the radius of gyration of a set of particles with (optional) radii and mass and (non-optional) coordinates. Either all particles must have mass or none of them.
Residue IMP::atom::get_residue | ( | Atom | d, |
bool | nothrow = false |
||
) |
Hierarchy IMP::atom::get_residue | ( | Hierarchy | mhd, |
unsigned int | index | ||
) |
Get the residue with the specified index.
Find the leaf containing the residue with the appropriate index. This is the PDB index, not the offset in the chain (if they are different).
The function returns a Hierarchy, rather than a Residue since the residue may not be explicitly represented and may just be part of some fragment.
ValueException | if mhd's type is not one of CHAIN, PROTEIN, NUCLEOTIDE |
ResidueType IMP::atom::get_residue_type | ( | char | c | ) |
Get the residue type from the 1-letter amino acid code.
ValueException | if an invalid character is passed. |
double IMP::atom::get_resolution | ( | kernel::Model * | m, |
kernel::ParticleIndex | pi | ||
) |
Return an estimate of the resolution of the hierarchy as used by Representation.
It is currently the inverse average radius of the leaves.
double IMP::atom::get_resolution | ( | Hierarchy | h | ) |
double IMP::atom::get_rigid_bodies_drms | ( | const Vector3DsOrXYZs0 & | m1, |
const Vector3DsOrXYZs1 & | m2, | ||
const IMP::IntRanges & | ranges | ||
) |
DRMS between to sets of rigid bodies. Points ij belonging to the same rigid body are not evaluated, as they are the same in both sets
[in] | m1 | set of points |
[in] | m2 | set of points |
[in] | ranges | of points considered to be the same rigid body. Eg, (0-29,30-49) means two rigid bodies, first with the 30 first atoms, second with the last 20 |
Definition at line 217 of file atom/distance.h.
double IMP::atom::get_rmsd | ( | const core::XYZs & | s0, |
const core::XYZs & | s1 | ||
) |
Get the rmsd between two lists of particles
double IMP::atom::get_rmsd | ( | const Vector3DsOrXYZs0 & | m1, |
const Vector3DsOrXYZs1 & | m2, | ||
const algebra::Transformation3D & | tr_for_second = algebra::get_identity_transformation_3d() |
||
) |
Get the rmsd between two lists of particles.
Definition at line 38 of file atom/distance.h.
double IMP::atom::get_rmsd | ( | const Selection & | s0, |
const Selection & | s1 | ||
) |
RMSD on a pair of Selections.
double IMP::atom::get_rmsd | ( | const Selection & | s0, |
const Selection & | s1, | ||
const algebra::Transformation3D & | tr_for_second | ||
) |
double IMP::atom::get_rmsd | ( | const algebra::Vector3Ds & | s0, |
const algebra::Vector3Ds & | s1, | ||
const algebra::Transformation3D & | tr_for_second | ||
) |
double IMP::atom::get_rmsd | ( | const core::XYZs & | s0, |
const core::XYZs & | s1, | ||
const IMP::algebra::Transformation3D & | tr_for_second | ||
) |
double IMP::atom::get_rmsd_transforming_first | ( | const IMP::algebra::Transformation3D & | tr, |
const core::XYZs & | s0, | ||
const core::XYZs & | s1 | ||
) |
Get the rmsd between two lists of particles
double IMP::atom::get_rmsd_transforming_first | ( | const algebra::Transformation3D & | tr, |
const Selection & | s0, | ||
const Selection & | s1 | ||
) |
RMSD on a pair of Selections.
Hierarchy IMP::atom::get_root | ( | Hierarchy | h | ) |
Return the root of the hierarchy.
Definition at line 425 of file atom/Hierarchy.h.
double IMP::atom::get_rotational_diffusion_coefficient | ( | const algebra::Rotation3Ds & | displacements, |
double | dt | ||
) |
Estimate the rotational diffusion coefficient of a particle from a list of displacements each taken after the given time step dt.
Float IMP::atom::get_secondary_structure_match_score | ( | SecondaryStructureResidue | ssr1, |
SecondaryStructureResidue | ssr2 | ||
) |
Compares the secondary structure probabilities of two residues.
int IMP::atom::get_state_index | ( | Hierarchy | h | ) |
Walk up the hierarchy to find the current state.
double IMP::atom::get_surface_area | ( | const Selection & | s | ) |
double IMP::atom::get_volume | ( | const Selection & | s | ) |
double IMP::atom::get_volume_from_mass | ( | double | m, |
ProteinDensityReference | ref = ALBER |
||
) |
Estimate the volume of a protein from its mass.
[in] | m | the mass for which we want to output the corresponding volume |
[in] | ref | the protein density reference used in the computation. As a default ref is the estimate published in Alber et. al, Structure 2005. |
double IMP::atom::get_volume_from_residue_type | ( | ResidueType | rt | ) |
Return an estimate for the volume of a given residue.
The volume estimates are taken from Pontius J, Richelle J, Wodak SJ., Deviations from standard atomic volumes as a quality measure for protein crystal structures, J Mol Biol. 1996 Nov 22;264(1):121-36.
ValueException | if a non-standard residue type is passed |
Hierarchy IMP::atom::read_mol2 | ( | base::TextInput | mol2_file, |
kernel::Model * | model, | ||
Mol2Selector * | mol2sel = nullptr |
||
) |
Create a hierarchy from a Mol2 file.
Hierarchies IMP::atom::read_multimodel_pdb | ( | base::TextInput | input, |
kernel::Model * | model, | ||
PDBSelector * | selector = get_default_pdb_selector() |
||
) |
Read all models from the PDB file.
Hierarchy IMP::atom::read_pdb | ( | base::TextInput | input, |
kernel::Model * | model, | ||
PDBSelector * | selector = get_default_pdb_selector() , |
||
bool | select_first_model = true |
||
) |
Read a all the molecules in the first model of the PDB file.
void IMP::atom::read_pdb | ( | base::TextInput | input, |
int | model, | ||
Hierarchy | h | ||
) |
Rewrite the coordinates of the passed hierarchy based on the contents of the first model in the PDB file.
The hierarchy must have been created by reading from a PDB file and the atom numbers must correspond between the files. These are not really checked.
A ValueException is thrown if there are insufficient models in the file.
core::RigidMember particles are handled by updating the core::RigidBody algebra::ReferenceFrame3D to align with the loaded particles. Bad things will happen if the loaded coordinates are not a rigid transform of the prior coordinates.
void IMP::atom::remove_charmm_untyped_atoms | ( | Hierarchy | hierarchy | ) |
Remove any atom from the Hierarchy that does not have a CHARMM type.
void IMP::atom::setup_as_approximation | ( | kernel::Particle * | h, |
const kernel::ParticlesTemp & | other | ||
) |
Set the mass, radius, residues, and coordinates to approximate the passed particles.
void IMP::atom::setup_as_approximation | ( | Hierarchy | h | ) |
Set the mass, radius, residues, and coordinates to approximate the passed particle based on the leaves of h.
SecondaryStructureResidue IMP::atom::setup_coarse_secondary_structure_residue | ( | const kernel::Particles & | ssr_ps, |
kernel::Model * | mdl, | ||
bool | winner_takes_all_per_res = false |
||
) |
Coarsen some SecondaryStructureResidues.
[in] | ssr_ps | The SSR-decorated particles to be combined |
[in] | mdl | The IMP kernel::Model |
[in] | winner_takes_all_per_res | Whether to set prob=1.0 for top scoring secondary structure type |
SecondaryStructureResidues IMP::atom::setup_coarse_secondary_structure_residues | ( | const kernel::Particles & | ssr_ps, |
kernel::Model * | mdl, | ||
int | coarse_factor, | ||
int | start_res_num, | ||
bool | winner_takes_all_per_res = false |
||
) |
Groups SecondaryStructureResidues into segments and then coarsens them. Useful if you have a long sequence and want to make several coarse nodes.
[in] | ssr_ps | The SSR-decorated particles to be combined |
[in] | mdl | The IMP kernel::Model |
[in] | coarse_factor | Group size |
[in] | start_res_num | Starting residue number for the provided sequence |
[in] | winner_takes_all_per_res | Whether to set prob=1.0 for top scoring secondary structure type |
void IMP::atom::show | ( | Hierarchy | h, |
std::ostream & | out = std::cout |
||
) |
void IMP::atom::show_molecular_hierarchy | ( | Hierarchy | h | ) |
void IMP::atom::transform | ( | atom::Hierarchy | h, |
const algebra::Transformation3D & | tr | ||
) |
Transform a hierarchy. This is aware of rigid bodies.
Transform a hierarchy, being aware of rigid bodies and intermediate nodes with coordinates. Rigid bodies must either be completely contained within the hierarchy or completely disjoint (no rigid bodies that contain particles in the hierarchy and other particles not in it are allowed).
void IMP::atom::write_mol2 | ( | Hierarchy | rhd, |
base::TextOutput | file_name | ||
) |
Write a ligand hierarchy as a mol2 file.
For now, this has to be a hierarchy created by read_mol2()
void IMP::atom::write_multimodel_pdb | ( | const Hierarchies & | mhd, |
base::TextOutput | out | ||
) |
Write the hierarchies one per frame.
void IMP::atom::write_pdb | ( | const Selection & | mhd, |
base::TextOutput | out, | ||
unsigned int | model = 1 |
||
) |
Write some atoms to a PDB.
void IMP::atom::write_pdb_of_c_alphas | ( | const Selection & | mhd, |
base::TextOutput | out, | ||
unsigned int | model = 1 |
||
) |
Write a hierarchy to a PDB as C_alpha atoms.
This method is used to write a non-atomic hierarchy into a PDB in a way that can be read by most programs. If the leaves are Residue particles then the index and residue type will be read from them. Otherwise default values will be used so that each leaf ends up in a separate residue.
const AtomType IMP::atom::AT_CA |
const ResidueType IMP::atom::GLY |
const ResidueType IMP::atom::UNK |
Unknown residue.