OpenStructure
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
Namespaces | Data Structures | Typedefs | Enumerations | Functions
ost::mol::alg Namespace Reference

Namespaces

namespace  helix_kinks
namespace  structure_analysis
namespace  superpose
namespace  trajectory_analysis
namespace  views

Data Structures

class  ClashingDistances
 List of reference atom-atom distances to detect clashes between non-bonded atoms. More...
class  StereoChemicalParams
 List of stereo chemical parameters (Bonds and angles) More...
class  UniqueAtomIdentifier
 Contains the infomation needed to uniquely identify an atom in a structure. More...
struct  SecStructureSegment
 a consecutive secondary structure element More...
struct  SuperpositionResult
 stores the number of cycles and the two final EntityViews (in case IterativeSuperposition was applied), the root mean square deviation of the two superposed objects and the transformation matrix (geom::Mat4) to superpose input1 onto input2. More...
class  SuperposerSVD
 efficiently superpose a bunch of models with the same number of atoms Choose either two EntityViews or two AtomViewLists. More...

Typedefs

typedef std::pair
< UniqueAtomIdentifier,
UniqueAtomIdentifier
UAtomIdentifiers
typedef std::map< std::pair
< UniqueAtomIdentifier,
UniqueAtomIdentifier >
, std::pair< Real, Real > > 
ResidueRDMap
typedef std::map
< ost::mol::ResNum,
ResidueRDMap
GlobalRDMap
typedef std::map
< UniqueAtomIdentifier, int > 
ExistenceMap
typedef std::vector
< SecStructureSegment
SecStructureSegments

Enumerations

enum  DensityType { HIGH_RESOLUTION, LOW_RESOLUTION }

Functions

Real DLLEXPORT_OST_MOL_ALG ClashScore (const mol::EntityView &ent_a, const mol::EntityView &ent_b)
Real DLLEXPORT_OST_MOL_ALG ClashScore (const mol::EntityHandle &ent_a, const mol::EntityView &ent_b)
Real DLLEXPORT_OST_MOL_ALG ClashScore (const mol::AtomHandle &atom, const mol::EntityView &ent_b)
Real DLLEXPORT_OST_MOL_ALG StericEnergy (const geom::Vec3 &pos1, Real r1, const geom::Vec3 &pos2, Real r2)
geom::Vec3 DLLEXPORT_OST_MOL_ALG CBetaPosition (const ResidueHandle &residue, Real bond_length=1.5)
geom::Vec3 DLLEXPORT_OST_MOL_ALG CBetaPosition (const geom::Vec3 &n_pos, const geom::Vec3 &ca_pos, const geom::Vec3 &c_pos, Real bond_length=1.5)
void DLLEXPORT_OST_MOL_ALG ConstructCBetas (EntityHandle &entity_handle, bool include_gly=false)
void DLLEXPORT_OST_MOL_ALG EntityToDensityScattering (const mol::EntityView &entity_view, img::MapHandle &map, Real falloff_start, Real falloff_end, bool clear_map_flag=false, Real source_wavelength=1.5418)
void DLLEXPORT_OST_MOL_ALG EntityToDensityRosetta (const mol::EntityView &entity_view, img::MapHandle &map, const DensityType &density_type, Real resolution, bool clear_map_flag=false, Real source_wavelength=1.5418)
ClashingDistances
DLLEXPORT_OST_MOL_ALG 
FillClashingDistances (std::vector< String > &stereo_chemical_props_file)
StereoChemicalParams
DLLEXPORT_OST_MOL_ALG 
FillStereoChemicalParams (const String &header, std::vector< String > &stereo_chemical_props_file)
EntityView DLLEXPORT_OST_MOL_ALG FilterClashes (const EntityView &ent, const ClashingDistances &min_distances, bool always_remove_bb=false)
EntityView DLLEXPORT_OST_MOL_ALG FilterClashes (const EntityHandle &ent, const ClashingDistances &min_distances, bool always_remove_bb=false)
EntityView DLLEXPORT_OST_MOL_ALG CheckStereoChemistry (const EntityView &ent, const StereoChemicalParams &bond_table, const StereoChemicalParams &angle_table, Real bond_tolerance, Real angle_tolerance, bool always_remove_bb=false)
EntityView DLLEXPORT_OST_MOL_ALG CheckStereoChemistry (const EntityHandle &ent, const StereoChemicalParams &bond_table, const StereoChemicalParams &angle_table, Real bond_tolerance, Real angle_tolerance, bool always_remove_bb=false)
std::pair< long int, long int >
DLLEXPORT_OST_MOL_ALG 
LocalDistDiffTest (const EntityView &mdl, const GlobalRDMap &dist_list, std::vector< Real > cutoff_list, int sequence_separation=0, const String &local_ldt_property_string="")
Real DLLEXPORT_OST_MOL_ALG LocalDistDiffTest (const EntityView &mdl, const EntityView &target, Real cutoff, Real max_dist, const String &local_ldt_property_string="")
Real DLLEXPORT_OST_MOL_ALG LocalDistDiffTest (const ost::seq::AlignmentHandle &aln, Real cutoff, Real max_dist, int ref_index=0, int mdl_index=1)
Real DLLEXPORT_OST_MOL_ALG LDDTHA (EntityView &v, const GlobalRDMap &global_dist_list, int sequence_separation=0)
GlobalRDMap CreateDistanceList (const EntityView &ref, Real max_dist)
GlobalRDMap CreateDistanceListFromMultipleReferences (const std::vector< EntityView > &ref_list, std::vector< Real > &cutoff_list, int sequence_separation, Real max_dist)
void PrintGlobalRDMap (const GlobalRDMap &glob_dist_list)
void PrintResidueRDMap (const ResidueRDMap &res_dist_list)
bool IsStandardResidue (String rn)
geom::Vec3List
DLLEXPORT_OST_MOL_ALG 
GetPosListFromView (const EntityView &view)
void DLLEXPORT_OST_MOL_ALG WrapEntityInPeriodicCell (EntityHandle eh, const geom::Vec3 cell_center, const geom::Vec3 basis_vec)
CoordGroupHandle
DLLEXPORT_OST_MOL_ALG 
SuperposeFrames (CoordGroupHandle &cg, EntityView &sel, int begin=0, int end=-1, int ref=-1)
CoordGroupHandle
DLLEXPORT_OST_MOL_ALG 
SuperposeFrames (CoordGroupHandle &cg, EntityView &sel, EntityView &ref_view, int begin=0, int end=-1)
SuperpositionResult
DLLEXPORT_OST_MOL_ALG 
SuperposeAtoms (const mol::AtomViewList &atoms1, const mol::AtomViewList &atoms2, bool apply_transform)
SuperpositionResult
DLLEXPORT_OST_MOL_ALG 
SuperposeSVD (const mol::EntityView &ev1, const mol::EntityView &ev2, bool apply_transform)
SuperpositionResult
DLLEXPORT_OST_MOL_ALG 
SuperposeSVD (const std::vector< geom::Vec3 > &pl1, const std::vector< geom::Vec3 > &pl2)
SuperpositionResult
DLLEXPORT_OST_MOL_ALG 
IterativeSuperposition (mol::EntityView &ev1, mol::EntityView &ev2, int ncycles, Real distance_threshold, bool apply_transform)
Real DLLEXPORT_OST_MOL_ALG CalculateRMSD (const mol::EntityView &ev1, const mol::EntityView &ev2, const geom::Mat4 &transformation)
geom::Vec3List
DLLEXPORT_OST_MOL_ALG 
AnalyzeAtomPos (const CoordGroupHandle &traj, const AtomHandle &a1, unsigned int stride=1)
geom::Vec3List
DLLEXPORT_OST_MOL_ALG 
AnalyzeCenterOfMassPos (const CoordGroupHandle &traj, const EntityView &sele, unsigned int stride=1)
std::vector< Real >
DLLEXPORT_OST_MOL_ALG 
AnalyzeDistanceBetwAtoms (const CoordGroupHandle &traj, const AtomHandle &a1, const AtomHandle &a2, unsigned int stride=1)
std::vector< Real >
DLLEXPORT_OST_MOL_ALG 
AnalyzeAngle (const CoordGroupHandle &traj, const AtomHandle &a1, const AtomHandle &a2, const AtomHandle &a3, unsigned int stride=1)
std::vector< Real >
DLLEXPORT_OST_MOL_ALG 
AnalyzeDistanceBetwCenterOfMass (const CoordGroupHandle &traj, const EntityView &sele1, const EntityView &sele2, unsigned int stride=1)
std::vector< Real >
DLLEXPORT_OST_MOL_ALG 
AnalyzeDihedralAngle (const CoordGroupHandle &traj, const AtomHandle &a1, const AtomHandle &a2, const AtomHandle &a3, const AtomHandle &a4, unsigned int stride=1)
std::vector< Real >
DLLEXPORT_OST_MOL_ALG 
AnalyzeRMSD (const CoordGroupHandle &traj, const EntityView &reference_view, const EntityView &sele, unsigned int stride=1)
Real DLLEXPORT_OST_MOL_ALG AnalyzeRMSF (const CoordGroupHandle &traj, const EntityView &selection, int from=0, int to=-1, unsigned int stride=1)
std::vector< Real >
DLLEXPORT_OST_MOL_ALG 
AnalyzeMinDistance (const CoordGroupHandle &traj, const EntityView &view1, const EntityView &view2, unsigned int stride=1)
std::vector< Real >
DLLEXPORT_OST_MOL_ALG 
AnalyzeMinDistanceBetwCenterOfMassAndView (const CoordGroupHandle &traj, const EntityView &view_cm, const EntityView &view_atoms, unsigned int stride=1)
std::vector< Real >
DLLEXPORT_OST_MOL_ALG 
AnalyzeAromaticRingInteraction (const CoordGroupHandle &traj, const EntityView &view_ring1, const EntityView &view_ring2, unsigned int stride=1)
void DLLEXPORT_OST_MOL_ALG AnalyzeAlphaHelixAxis (const CoordGroupHandle &traj, const EntityView &prot_seg, geom::Vec3List &directions, geom::Vec3List &centers, unsigned int stride=1)
void DLLEXPORT_OST_MOL_ALG AnalyzeBestFitLine (const CoordGroupHandle &traj, const EntityView &prot_seg, geom::Vec3List &directions, geom::Vec3List &centers, unsigned int stride=1)
void DLLEXPORT_OST_MOL_ALG AnalyzeBestFitPlane (const CoordGroupHandle &traj, const EntityView &prot_seg, geom::Vec3List &normals, geom::Vec3List &origins, unsigned int stride=1)
EntityHandle DLLEXPORT_OST_MOL_ALG CreateMeanStructure (const CoordGroupHandle &traj, const EntityView &selection, int from=0, int to=-1, unsigned int stride=1)
std::vector< Real >
DLLEXPORT_OST_MOL_ALG 
AnalyzeHelicity (const CoordGroupHandle &traj, const EntityView &prot_seg, unsigned int stride=1)
def FillClashingDistancesFromFile
def FillBondStereoChemicalParamsFromFile
def FillAngleStereoChemicalParamsFromFile
def DefaultClashingDistances
def DefaultBondStereoChemicalParams
def DefaultAngleStereoChemicalParams

Typedef Documentation

typedef std::map<UniqueAtomIdentifier,int> ExistenceMap

Definition at line 95 of file local_dist_diff_test.hh.

Global distance list.

Container for all the residue-based interatomic distance lists that are checked in a Local Distance Difference Test and belong to the same structure

Definition at line 92 of file local_dist_diff_test.hh.

typedef std::map<std::pair<UniqueAtomIdentifier,UniqueAtomIdentifier>,std::pair<Real,Real> > ResidueRDMap

Residue distance list.

Container for all the interatomic distances that are checked in a Local Distance Difference Test and are originating from a single specific residue

Definition at line 86 of file local_dist_diff_test.hh.

Definition at line 46 of file sec_structure_segments.hh.

Definition at line 80 of file local_dist_diff_test.hh.


Enumeration Type Documentation

type of density being created by the EntityToDensity function

Enumerator:
HIGH_RESOLUTION 
LOW_RESOLUTION 

Definition at line 13 of file entity_to_density.hh.


Function Documentation

void DLLEXPORT_OST_MOL_ALG ost::mol::alg::AnalyzeAlphaHelixAxis ( const CoordGroupHandle &  traj,
const EntityView &  prot_seg,
geom::Vec3List directions,
geom::Vec3List centers,
unsigned int  stride = 1 
)
std::vector<Real> DLLEXPORT_OST_MOL_ALG ost::mol::alg::AnalyzeAngle ( const CoordGroupHandle &  traj,
const AtomHandle &  a1,
const AtomHandle &  a2,
const AtomHandle &  a3,
unsigned int  stride = 1 
)
std::vector<Real> DLLEXPORT_OST_MOL_ALG ost::mol::alg::AnalyzeAromaticRingInteraction ( const CoordGroupHandle &  traj,
const EntityView &  view_ring1,
const EntityView &  view_ring2,
unsigned int  stride = 1 
)
geom::Vec3List DLLEXPORT_OST_MOL_ALG ost::mol::alg::AnalyzeAtomPos ( const CoordGroupHandle &  traj,
const AtomHandle &  a1,
unsigned int  stride = 1 
)
void DLLEXPORT_OST_MOL_ALG ost::mol::alg::AnalyzeBestFitLine ( const CoordGroupHandle &  traj,
const EntityView &  prot_seg,
geom::Vec3List directions,
geom::Vec3List centers,
unsigned int  stride = 1 
)
void DLLEXPORT_OST_MOL_ALG ost::mol::alg::AnalyzeBestFitPlane ( const CoordGroupHandle &  traj,
const EntityView &  prot_seg,
geom::Vec3List normals,
geom::Vec3List origins,
unsigned int  stride = 1 
)
geom::Vec3List DLLEXPORT_OST_MOL_ALG ost::mol::alg::AnalyzeCenterOfMassPos ( const CoordGroupHandle &  traj,
const EntityView &  sele,
unsigned int  stride = 1 
)
std::vector<Real> DLLEXPORT_OST_MOL_ALG ost::mol::alg::AnalyzeDihedralAngle ( const CoordGroupHandle &  traj,
const AtomHandle &  a1,
const AtomHandle &  a2,
const AtomHandle &  a3,
const AtomHandle &  a4,
unsigned int  stride = 1 
)
std::vector<Real> DLLEXPORT_OST_MOL_ALG ost::mol::alg::AnalyzeDistanceBetwAtoms ( const CoordGroupHandle &  traj,
const AtomHandle &  a1,
const AtomHandle &  a2,
unsigned int  stride = 1 
)
std::vector<Real> DLLEXPORT_OST_MOL_ALG ost::mol::alg::AnalyzeDistanceBetwCenterOfMass ( const CoordGroupHandle &  traj,
const EntityView &  sele1,
const EntityView &  sele2,
unsigned int  stride = 1 
)
std::vector<Real> DLLEXPORT_OST_MOL_ALG ost::mol::alg::AnalyzeHelicity ( const CoordGroupHandle &  traj,
const EntityView &  prot_seg,
unsigned int  stride = 1 
)
std::vector<Real> DLLEXPORT_OST_MOL_ALG ost::mol::alg::AnalyzeMinDistance ( const CoordGroupHandle &  traj,
const EntityView &  view1,
const EntityView &  view2,
unsigned int  stride = 1 
)
std::vector<Real> DLLEXPORT_OST_MOL_ALG ost::mol::alg::AnalyzeMinDistanceBetwCenterOfMassAndView ( const CoordGroupHandle &  traj,
const EntityView &  view_cm,
const EntityView &  view_atoms,
unsigned int  stride = 1 
)
std::vector<Real> DLLEXPORT_OST_MOL_ALG ost::mol::alg::AnalyzeRMSD ( const CoordGroupHandle &  traj,
const EntityView &  reference_view,
const EntityView &  sele,
unsigned int  stride = 1 
)
Real DLLEXPORT_OST_MOL_ALG ost::mol::alg::AnalyzeRMSF ( const CoordGroupHandle &  traj,
const EntityView &  selection,
int  from = 0,
int  to = -1,
unsigned int  stride = 1 
)
Real DLLEXPORT_OST_MOL_ALG ost::mol::alg::CalculateRMSD ( const mol::EntityView &  ev1,
const mol::EntityView &  ev2,
const geom::Mat4 transformation 
)

calculates RMSD for two entity view

geom::Vec3 DLLEXPORT_OST_MOL_ALG ost::mol::alg::CBetaPosition ( const ResidueHandle &  residue,
Real  bond_length = 1.5 
)
geom::Vec3 DLLEXPORT_OST_MOL_ALG ost::mol::alg::CBetaPosition ( const geom::Vec3 n_pos,
const geom::Vec3 ca_pos,
const geom::Vec3 c_pos,
Real  bond_length = 1.5 
)
EntityView DLLEXPORT_OST_MOL_ALG ost::mol::alg::CheckStereoChemistry ( const EntityView &  ent,
const StereoChemicalParams &  bond_table,
const StereoChemicalParams &  angle_table,
Real  bond_tolerance,
Real  angle_tolerance,
bool  always_remove_bb = false 
)

Filters a structure based on detected stereo-chemical violations. Entity version.

If a stereo-chemical violation (i.e., a bond or an angle with a value outside the range defined by the mean value, the standard deviation and the tolerance parameter) is detected in a residue's side-chain, all atoms in the side chain are removed from the structure. If a violation is detected in the backbone, all atoms in the residue are removed. This behavior is changed by the always_remove_bb flag: when the flag is set to true all atoms in the residue are removed even if a violation is just detected in the side-chain

EntityView DLLEXPORT_OST_MOL_ALG ost::mol::alg::CheckStereoChemistry ( const EntityHandle &  ent,
const StereoChemicalParams &  bond_table,
const StereoChemicalParams &  angle_table,
Real  bond_tolerance,
Real  angle_tolerance,
bool  always_remove_bb = false 
)

Filters a structure based on detected stereo-chemical violations. Handle version.

If a stereo-chemical violation (i.e., a bond or an angle with a value outside the range defined by the mean value, the standard deviation and the tolerance parameter) is detected in a residue's side-chain, all atoms in the side chain are removed from the structure. If a violation is detected in the backbone, all atoms in the residue are removed. This behavior is changed by the always_remove_bb flag: when the flag is set to true all atoms in the residue are removed even if a violation is just detected in the side-chain

void DLLEXPORT_OST_MOL_ALG ost::mol::alg::ConstructCBetas ( EntityHandle &  entity_handle,
bool  include_gly = false 
)
GlobalRDMap ost::mol::alg::CreateDistanceList ( const EntityView &  ref,
Real  max_dist 
)

Creates a list of distances to check during a Local Difference Distance Test.

Requires a reference structure and an inclusion radius (max_dist)

GlobalRDMap ost::mol::alg::CreateDistanceListFromMultipleReferences ( const std::vector< EntityView > &  ref_list,
std::vector< Real > &  cutoff_list,
int  sequence_separation,
Real  max_dist 
)

Creates a list of distances to check during a Local Difference Distance Test starting from multiple reference structures.

Requires a list of reference structure and an inclusion radius (max_dist).

The structures in the list have to be properly prepared before being passed to the function. Corresponding residues in the structures must have the same residue number, the same chain name, etc. Gaps are allowed and automatically dealt with: if information about a distance is present in at least one of the structures, it will be considered.

If a distance between two atoms is shorter than the inclusion radius in all structures in which the two atoms are present, it is included in the list. However, if the distance is longer than the inclusion radius in at least one of the structures, it is not be considered a local interaction and is exluded from the list

The function takes care of residues with ambigous symmetric sidechains. To decide which naming convention to use, the functions computes a local distance score of each reference structure with the first reference structure in the list, using only non ambigously-named atoms. It picks then the naming convention that gives the highest score, guaranteeing that all references are processed with the correct atom names.

The cutoff list that will later be used to compute the Local Distance Difference Test score and the sequence separation parameter must be passed to the function. These parameters do not influence the output distance list, which always includes all distances within the provided max_dist (to make it consistent with the single-reference corresponding function). However, the parameters are used when dealing with the naming convention of residues with ambiguous nomenclature.

EntityHandle DLLEXPORT_OST_MOL_ALG ost::mol::alg::CreateMeanStructure ( const CoordGroupHandle &  traj,
const EntityView &  selection,
int  from = 0,
int  to = -1,
unsigned int  stride = 1 
)
def ost.mol.alg.DefaultAngleStereoChemicalParams ( )

Definition at line 48 of file __init__.py.

def ost.mol.alg.DefaultBondStereoChemicalParams ( )

Definition at line 39 of file __init__.py.

def ost.mol.alg.DefaultClashingDistances ( )

Definition at line 30 of file __init__.py.

void DLLEXPORT_OST_MOL_ALG ost::mol::alg::EntityToDensityRosetta ( const mol::EntityView &  entity_view,
img::MapHandle &  map,
const DensityType &  density_type,
Real  resolution,
bool  clear_map_flag = false,
Real  source_wavelength = 1.5418 
)

create a density representation of an entity in a density map

This function creates a density representation of the entity provided by the user in a density map, also provided by the user. The user can choose the type of density of the output map:

ROSETTA_HIGH_RESOLUTION gaussian spheres in real space to represent density, one per atom, see Dimaio et al., Refinement of Protein Structures into Low-Resolution Density Maps Using Rosetta. Journal of Molecular Biology (2009) pp. 1-10 ROSETTA_LOW_RESOLUTION guassian spheres in real space to represent density, one per residue. See reference above. Only useful at low resolution

The user can also choose if the density map should be cleared of its previous content before creating the density representation.

The user must also provide a resolution parameter.

This function will only create a density represenation of the entities (or portion of entities ) that fall within the borders of the map. The user must take care that this condition is verified for all entities for which he wants a representation.

void DLLEXPORT_OST_MOL_ALG ost::mol::alg::EntityToDensityScattering ( const mol::EntityView &  entity_view,
img::MapHandle &  map,
Real  falloff_start,
Real  falloff_end,
bool  clear_map_flag = false,
Real  source_wavelength = 1.5418 
)

create a density representation of an entity in a density map (using electron scattering factors)

This functions creates a density representation of the entity provided by the user in a density map, in Fourier space using the correct scattering factors for the elements involved.

The user can also choose if the density map should be cleared of its previous content before creating the density representation.

The density is generated in Fourier space. In order to avoid artifacts in the final density representation, the function avoids sharp frequency cutoffs by applying a Gaussian falloff. The user must provide the resolutions at which the cutoff should begin and end, as opposed to a single resolution cutoff value.

This function will only create a density represenation of the entities (or portion of entities ) that fall within the borders of the map. The user must take care that this condition is verified for all entities for which he wants a representation.

def ost.mol.alg.FillAngleStereoChemicalParamsFromFile (   filename)

Definition at line 23 of file __init__.py.

def ost.mol.alg.FillBondStereoChemicalParamsFromFile (   filename)

Definition at line 16 of file __init__.py.

ClashingDistances DLLEXPORT_OST_MOL_ALG ost::mol::alg::FillClashingDistances ( std::vector< String > &  stereo_chemical_props_file)

Fills a list of reference clashing distances from the content of a parameter file.

Requires a list of strings holding the contents of a parameter file, one line per string

def ost.mol.alg.FillClashingDistancesFromFile (   filename,
  default_clashing_distance,
  default_tolerance 
)

Definition at line 9 of file __init__.py.

StereoChemicalParams DLLEXPORT_OST_MOL_ALG ost::mol::alg::FillStereoChemicalParams ( const String header,
std::vector< String > &  stereo_chemical_props_file 
)

Fills a list of stereo-chemical statistics from the content of a parameter file.

Requires a list of strings holding the contents of a parameter file, one line per string The header can be 'Bonds' to read bond statistics or 'Angles' to read angle statistics

EntityView DLLEXPORT_OST_MOL_ALG ost::mol::alg::FilterClashes ( const EntityView &  ent,
const ClashingDistances &  min_distances,
bool  always_remove_bb = false 
)

Filters a structure based on detected clashes between non bonded atoms. Entity version.

If a clash between two atoms (distance shorter than reference clashing distance - tolerance threshold) is detected in a residue's side-chain, all atoms in the side chain are removed from the structure If a clash is detected in the backbone, all atoms in the residue are removed. This behavior is changed by the always_remove_bb flag: when the flag is set to true all atoms in the residue are removed even if a clash is just detected in the side-chain

EntityView DLLEXPORT_OST_MOL_ALG ost::mol::alg::FilterClashes ( const EntityHandle &  ent,
const ClashingDistances &  min_distances,
bool  always_remove_bb = false 
)

Filters a structure based on detected clashes between non bonded atoms. Handle version.

If a clash between two atoms (distance shorter than reference clashing distance - tolerance threshold) is detected in a residue's side-chain, all atoms in the side chain are removed from the structure If a clash is detected in the backbone, all atoms in the residue are removed. This behavior is changed by the always_remove_bb flag: when the flag is set to true all atoms in the residue are removed even if a clash is just detected in the side-chain

geom::Vec3List DLLEXPORT_OST_MOL_ALG ost::mol::alg::GetPosListFromView ( const EntityView &  view)
bool ost::mol::alg::IsStandardResidue ( String  rn)
SuperpositionResult DLLEXPORT_OST_MOL_ALG ost::mol::alg::IterativeSuperposition ( mol::EntityView &  ev1,
mol::EntityView &  ev2,
int  ncycles,
Real  distance_threshold,
bool  apply_transform 
)

iterative superposition

Real DLLEXPORT_OST_MOL_ALG ost::mol::alg::LDDTHA ( EntityView &  v,
const GlobalRDMap &  global_dist_list,
int  sequence_separation = 0 
)

Computes the Local Distance Difference High-Accuracy Test given a list of distances to check.

Computes the Local Distance Difference High-Accuracy Test (with threshold 0.5,1,2 and 4 Angstrom) Requires a list of distances to check and a model for which the score is computed

A sequence separation parameter can be passed to the function. If this happens, only distances between residues whose separation is higher than the provided parameter are considered when computing the score.

std::pair<long int,long int> DLLEXPORT_OST_MOL_ALG ost::mol::alg::LocalDistDiffTest ( const EntityView &  mdl,
const GlobalRDMap &  dist_list,
std::vector< Real cutoff_list,
int  sequence_separation = 0,
const String local_ldt_property_string = "" 
)

Calculates number of distances conserved in a model, given a list of distances to check and a model.

Calculates the two values needed to determine the Local Distance Difference Test for a given model, i.e. the number of conserved distances in the model and the number of total distances in the reference structure. The function requires a list of distances to check, a model on which the distances are checked, and a list of tolerance thresholds that are used to determine if the distances are conserved.

The function only processes standard residues in the first chains of the model and of the reference For residues with symmetric sidechains (GLU, ASP, ARG, VAL, PHE, TYR), the naming of the atoms is ambigous. For these residues, the function computes the Local Distance Difference Test score that each naming convention would generate when considering all non-ambigous surrounding atoms. The solution that gives higher score is then picked to compute the final Local Difference Distance Score for the whole model.

A sequence separation parameter can be passed to the function. If this happens, only distances between residues whose separation is higher than the provided parameter are considered when computing the score.

If a string is provided as an argument to the function, residue-per-residue statistics are stored as residue properties. Specifically, the local residue-based lddt score is stored in a float property named as the provided string, while the residue-based number of conserved and total distances are saved in two int properties named <string>_conserved and <string>_total.

Real DLLEXPORT_OST_MOL_ALG ost::mol::alg::LocalDistDiffTest ( const EntityView &  mdl,
const EntityView &  target,
Real  cutoff,
Real  max_dist,
const String local_ldt_property_string = "" 
)

Calculates the Local Distance Difference Score for a given model with respect to a given target.

Calculates the Local Distance Difference Test score for a given model with respect to a given reference structure. Requires a model, a reference structure, a list of thresholds that are used to determine if distances are conserved, and an inclusion radius value used to determine which distances are checked.

The function only processes standard residues in the first chains of the model and of the reference For residues with symmetric sidechains (GLU, ASP, ARG, VAL, PHE, TYR), the naming of the atoms is ambigous. For these residues, the function computes the Local Distance Difference Test score that each naming convention would generate when considering all non-ambigous surrounding atoms. The solution that gives higher score is then picked to compute the final Local Difference Distance Score for the whole model.

If a string is provided as an argument to the function, residue-per-residue statistics are stored as residue properties. Specifically, the local residue-based lddt score is stored in a float property named as the provided string, while the residue-based number of conserved and total distances are saved in two int properties named <string>_conserved and <string>_total.

Real DLLEXPORT_OST_MOL_ALG ost::mol::alg::LocalDistDiffTest ( const ost::seq::AlignmentHandle aln,
Real  cutoff,
Real  max_dist,
int  ref_index = 0,
int  mdl_index = 1 
)

Calculates the Local Distance Difference Test score for a given model starting from an alignment between a reference structure and the model.

Calculates the Local Distance Difference Test score given an alignment between a model and a taget structure. Requires a threshold on which to calculate the score and an inclusion radius to determine the interatiomic distances to check. Obviously, the strucvtures of the model and the reference must be attached to the alignment. By default the first structure in the alignment is considered the reference and the second is considered the model, but this can be changed by passing to the function the indexes of the two structures in the ref_index and mdl_index parameters. BEWARE: This algorithm uses the old version of the Local Distance Difference Test (multiple cycles, single threshold, etc. ) and will give a slightly different result than the other functions

void ost::mol::alg::PrintGlobalRDMap ( const GlobalRDMap &  glob_dist_list)

Prints all distances in a global distance list to standard output.

void ost::mol::alg::PrintResidueRDMap ( const ResidueRDMap &  res_dist_list)

Prints all distances in a residue distance list to standard output.

SuperpositionResult DLLEXPORT_OST_MOL_ALG ost::mol::alg::SuperposeAtoms ( const mol::AtomViewList &  atoms1,
const mol::AtomViewList &  atoms2,
bool  apply_transform 
)

takes the corresponding atoms and superposes them

CoordGroupHandle DLLEXPORT_OST_MOL_ALG ost::mol::alg::SuperposeFrames ( CoordGroupHandle &  cg,
EntityView &  sel,
int  begin = 0,
int  end = -1,
int  ref = -1 
)

returns a superposed version of coord group, superposed on a reference frame

CoordGroupHandle DLLEXPORT_OST_MOL_ALG ost::mol::alg::SuperposeFrames ( CoordGroupHandle &  cg,
EntityView &  sel,
EntityView &  ref_view,
int  begin = 0,
int  end = -1 
)

returns a superposed version of coord group, superposed on a reference view

SuperpositionResult DLLEXPORT_OST_MOL_ALG ost::mol::alg::SuperposeSVD ( const mol::EntityView &  ev1,
const mol::EntityView &  ev2,
bool  apply_transform 
)

superposes two entity views

SuperpositionResult DLLEXPORT_OST_MOL_ALG ost::mol::alg::SuperposeSVD ( const std::vector< geom::Vec3 > &  pl1,
const std::vector< geom::Vec3 > &  pl2 
)

superposes two pointlists

void DLLEXPORT_OST_MOL_ALG ost::mol::alg::WrapEntityInPeriodicCell ( EntityHandle  eh,
const geom::Vec3  cell_center,
const geom::Vec3  basis_vec 
)