ProSHADE
0.6.5 (NOV 2018)
Protein Shape Descriptors and Symmetry Detection
|
This is the executive class responsible for comparing two or more structures. More...
#include <ProSHADE_internal.h>
Public Member Functions | |
ProSHADE_compareOneAgainstAll (ProSHADE_data *oneStr, std::vector< ProSHADE_data *> *allStrs, std::vector< int > ignoreL, double matrixPowerWeight, int verbose) | |
Contructor for the ProSHADE_compareOneAgainstAll class. More... | |
ProSHADE_compareOneAgainstAll (ProSHADE_data *oneStr, ProSHADE_data *twoStr, std::vector< ProSHADE_data *> *allStrs, std::vector< int > ignoreL, double matrixPowerWeight, int verbose) | |
Contructor for the ProSHADE_compareOneAgainstAll class with reverse phase allowed. More... | |
~ProSHADE_compareOneAgainstAll () | |
Destructor for the ProSHADE_compareOneAgainstAll class. More... | |
std::vector< double > | getEnergyLevelsDistance (int verbose, ProSHADE::ProSHADE_settings *settings) |
This function computes the energy level descriptor value from the first structure to all remaining structures. More... | |
void | precomputeTrSigmaDescriptor (double shellSpacing, std::vector< unsigned int > *glIntegOrderVec, ProSHADE::ProSHADE_settings *settings) |
This function computes the E matrices required for the trace sigma descriptor, the rotation function descriptor and SO3 transform. More... | |
std::vector< double > | getTrSigmaDistance (int verbose, ProSHADE::ProSHADE_settings *settings) |
This function computes the trace sigma descriptor distances. More... | |
void | getSO3InverseMap (ProSHADE::ProSHADE_settings *settings) |
This function is responsible for computing the SO3 inverse transform. More... | |
std::vector< std::array< double, 3 > > | getEulerAngles (ProSHADE::ProSHADE_settings *settings) |
This function finds the highest peak in the SO3 inverse transform map and sets it as the optimal overlay angle. More... | |
void | generateWignerMatrices (ProSHADE::ProSHADE_settings *settings) |
This function is responsible for computing the Wigner D matrices for full rotation function distance computation purposes. More... | |
std::vector< double > | getRotCoeffDistance (int verbose, ProSHADE::ProSHADE_settings *settings) |
This function computes the full rotation function descriptor distances. More... | |
void | alignDensities (ProSHADE::ProSHADE_settings *settings) |
Takes the internal objects with and without phases and aligns them to all the other objects, which should be read with reverse phased database. More... | |
Friends | |
class | ProSHADE_distances |
This class is used to compare two or more structures and is currently used for distances computation, but should not really be used for directly by the user, as there is designated class for this (ProSHADE_distances) and this class uses the current and best tested approach.
Definition at line 472 of file ProSHADE_internal.h.
ProSHADE_internal::ProSHADE_compareOneAgainstAll::ProSHADE_compareOneAgainstAll | ( | ProSHADE_data * | oneStr, |
std::vector< ProSHADE_data *> * | allStrs, | ||
std::vector< int > | ignoreL, | ||
double | matrixPowerWeight, | ||
int | verbose | ||
) |
This is the constructor for the executive class responsible for comparing two or more structures with different settings. This is currently internally used for distances computation; however the user should use the designated distances computation class ProSHADE_distances. The constructor is responsible for sanity checks, setting up parameters and determining the integration variables.
[in] | one | This is the pointer to the first structure data holding object which will be compared againt the rest. |
[in] | all | This is the pointer to an array of structure data containing all the other structures to be compared against the one. |
[in] | ignoreL | This vector should contain all bands to be ignored in the computations. |
[in] | matrixPowerWeight | This is the weighting parameter for energy level distances computation. |
[out] | X | Object capable of comparing the two or more structures. |
Definition at line 801 of file ProSHADE_internal.cpp.
ProSHADE_internal::ProSHADE_compareOneAgainstAll::ProSHADE_compareOneAgainstAll | ( | ProSHADE_data * | oneStr, |
ProSHADE_data * | twoStr, | ||
std::vector< ProSHADE_data *> * | allStrs, | ||
std::vector< int > | ignoreL, | ||
double | matrixPowerWeight, | ||
int | verbose | ||
) |
This is the constructor for the executive class responsible for comparing two or more structures with different settings. This is currently internally used for distances computation; however the user should use the designated distances computation class ProSHADE_distances. The constructor is responsible for sanity checks, setting up parameters and determining the integration variables.
[in] | one | This is the pointer to the first structure data holding object which will be compared againt the rest. |
[in] | one | This is the pointer to the second structure data holding object (same as the first, but with reverse phase) which will be compared againt the rest. |
[in] | all | This is the pointer to an array of structure data containing all the other structures to be compared against the one. This assumes this vector has been read from a database with phase reversal allowed, so that one structure is with phase and the other is without. Othewise, there be dragons. |
[in] | ignoreL | This vector should contain all bands to be ignored in the computations. |
[in] | matrixPowerWeight | This is the weighting parameter for energy level distances computation. |
[out] | X | Object capable of comparing the two or more structures. |
Definition at line 865 of file ProSHADE_internal.cpp.
ProSHADE_internal::ProSHADE_compareOneAgainstAll::~ProSHADE_compareOneAgainstAll | ( | ) |
This destructor is empty as its role has been taken over by the classes this are directly used by the user; this is because the data should only be deleted when the user is done using the user-accessible class and not before that.
Definition at line 932 of file ProSHADE_internal.cpp.
void ProSHADE_internal::ProSHADE_compareOneAgainstAll::alignDensities | ( | ProSHADE::ProSHADE_settings * | settings | ) |
This function is called if a structure is compared to a database with reverse phased (i.e. with both phased and phaseless data) structures. It will use both of the data types to find the optimal translation and rotation and then it will apply these to the structures, releasing the old intertwined vector and returning a new vector of aligned structures. This is a very particular functions, use only for the purpose specified.
[in] | settings | A settings object having all the required values. |
Definition at line 6085 of file ProSHADE_SO3transform.cpp.
void ProSHADE_internal::ProSHADE_compareOneAgainstAll::generateWignerMatrices | ( | ProSHADE::ProSHADE_settings * | settings | ) |
This function makes use of the SOFT2.0 package functions and logic to compute the Wigner d matrices and consequently the Wigner D matrices, which it stores internally in the distances computation class. The main purpose of the Wigner D matrices is to rotate the spherical harmonics coefficients.
[in] | settings | The ProSHADE_settings container class with information about how to write the HTML report. |
Definition at line 207 of file ProSHADE_wigner.cpp.
std::vector< double > ProSHADE_internal::ProSHADE_compareOneAgainstAll::getEnergyLevelsDistance | ( | int | verbose, |
ProSHADE::ProSHADE_settings * | settings | ||
) |
This function computes the Pearson's correlation coefficient between the weighted RRP matrices of two structures and proceeds to compute the energy level descriptor from these.
[in] | verbose | Int value stating how loud standard output should be. |
[in] | settings | The settings contained with inforamation about printing into the HTML report file, if need be. |
[out] | X | A vector of energy level distances from the first structure to all the other structeres. |
Definition at line 310 of file ProSHADE_descriptors.cpp.
std::vector< std::array< double, 3 > > ProSHADE_internal::ProSHADE_compareOneAgainstAll::getEulerAngles | ( | ProSHADE::ProSHADE_settings * | settings | ) |
This function simply finds the highest peak in the SO3 inverse transform map, converst the coordinates to ZXZ Euler angles and sets these as the optimal overlay angles. It does not check for the zero-angle peaks, as it assumes the two strucrures are different in which case the zero-angle peaks are valied as well. It should be used only for distances computation and never for symmetry detection purposes.
[in] | settings | The settings class container with information about how to print the HTML report. |
[out] | X | This function returns a vector of Euler angles (array of three numbers) for all structures which are compared to the single structure in multiple structure comparison. |
Definition at line 516 of file ProSHADE_SO3transform.cpp.
std::vector< double > ProSHADE_internal::ProSHADE_compareOneAgainstAll::getRotCoeffDistance | ( | int | verbose, |
ProSHADE::ProSHADE_settings * | settings | ||
) |
Since it would require one extra array the size of the spherical harmonics array to hold the Wigner D * C (C is spherical harmonics coefficients matrix) matrix multiplication result, this function takes the already available E matrices and multiplies these with the Wigner D matrices, this saving space. Otherwise, this is the rotation function known from molecular replacement and the weighted sum of the results is then the distance.
[in] | verbose | Int value stating how loud standard output should be. |
[in] | settings | The ProSHADE_settings class instance specifying how to write the HTML report output. |
[out] | X | A vector of double precision numbers which are the distances from the first structure to all remaining structures. |
Definition at line 1063 of file ProSHADE_descriptors.cpp.
void ProSHADE_internal::ProSHADE_compareOneAgainstAll::getSO3InverseMap | ( | ProSHADE::ProSHADE_settings * | settings | ) |
The SO3 transform and its inverse are required by the distances computing class to obtain the best angles at which the two structures overlay, so that the full rotation function distance could be computed using these.
[in] | settings | The settings class container with information about how to print the HTML report. |
Definition at line 217 of file ProSHADE_SO3transform.cpp.
std::vector< double > ProSHADE_internal::ProSHADE_compareOneAgainstAll::getTrSigmaDistance | ( | int | verbose, |
ProSHADE::ProSHADE_settings * | settings | ||
) |
This function basically takes the E matrix and calls LAPACK to decompose it onto the singular values. The returned singular values matrix then has its diagonal summed and this is the trace sigma descriptor in simplicity. The function does several other minor things, but they are not of much interest here.
[in] | verbose | Int value stating how loud standard output should be. |
[in] | settings | ProSHADE_settings container class with information about the HTML report printing, if needed. |
[out] | X | A vector of double precision numbers which are the distances from the first to all remaining structures. |
Definition at line 879 of file ProSHADE_descriptors.cpp.
void ProSHADE_internal::ProSHADE_compareOneAgainstAll::precomputeTrSigmaDescriptor | ( | double | shellSpacing, |
std::vector< unsigned int > * | glIntegOrderVec, | ||
ProSHADE::ProSHADE_settings * | settings | ||
) |
This function computes the integral over r for the c1_{l,m} * c2_{l,m} spherical harmonics values. These are used by all three of the following: the trace sigma descriptor, the SO3 transform and the full rotation function. As such, there are two copies of this function, one for the symmetry detection class and one for the distances computation; this is the version for the distances computation.
[in] | shellSpacing | How far apart individual shells are in angstroms? |
[in] | glIntegOrderVec | Pointer to the vector of integration values for multiple structure comparisons. |
[in] | settings | ProSHADE_settings container with infotmation about HTML reporting. |
Definition at line 646 of file ProSHADE_descriptors.cpp.