ProSHADE
0.6.5 (NOV 2018)
Protein Shape Descriptors and Symmetry Detection
|
This is the class which computes the symmetry in a single structure. More...
#include <ProSHADE_internal.h>
Public Member Functions | |
ProSHADE_symmetry (ProSHADE::ProSHADE_settings *settings) | |
Contructor for the ProSHADE_symmetry class. More... | |
ProSHADE_symmetry () | |
Contructor for the ProSHADE_symmetry class. More... | |
void | printResultsClear (int verbose) |
This function prints the cleared results to the screen. More... | |
void | printResultsClearHTML (ProSHADE::ProSHADE_settings *settings) |
This function prints the cleared results to the HTML file. More... | |
void | printResultsRequest (std::string symmetryType, unsigned int symmetryFold, int verbose) |
This function prints the cleared results to the screen. More... | |
void | printResultsRequestHTML (std::string symmetryType, unsigned int symmetryFold, ProSHADE::ProSHADE_settings *settings) |
This function prints the cleared results to the HTML report file. More... | |
std::vector< std::array< double, 8 > > | getRotFnPeaks () |
This function gives the user programmatical access to the symmetry peaks of the symmetry detection procedure. More... | |
std::vector< std::array< double, 5 > > | getCSymmetries () |
This function gives the user programmatical access to the detected C-symmetries. More... | |
std::vector< std::vector< std::array< double, 6 > > > | getDSymmetries () |
This function gives the user programmatical access to the detected D-symmetries. More... | |
std::vector< std::array< double, 5 > > | generateTetrAxes (ProSHADE_comparePairwise *cmpObj, ProSHADE::ProSHADE_settings *settings, std::vector< std::array< double, 5 > > tetrSymm, std::vector< std::array< double, 5 > > allCs, double axisErrorTolerance=0.1, int verbose=0) |
This function generates all the tetrahedral symmetry axes from a pair detected by findTetrSymmetry function. More... | |
std::vector< std::array< double, 5 > > | generateTetrElements (std::vector< std::array< double, 5 > > symmAxes, ProSHADE::ProSHADE_settings *settings, int verbose=0) |
This function generates the 12 unique tetrahedral symmetry group elements from the already detected axes. More... | |
std::vector< std::array< double, 5 > > | generateOctaAxes (ProSHADE_comparePairwise *cmpObj, ProSHADE::ProSHADE_settings *settings, std::vector< std::array< double, 5 > > octaSymm, std::vector< std::array< double, 5 > > allCs, double axisErrorTolerance=0.1, int verbose=0) |
This function generates all the ctahedral symmetry elements from a pair detected by findOctaSymmetry function. More... | |
std::vector< std::array< double, 5 > > | generateOctaElements (std::vector< std::array< double, 5 > > symmAxes, ProSHADE::ProSHADE_settings *settings, int verbose=0) |
This function generates the 24 octahedral symmetry group elements from the axes generated by generateOctaAxes function. More... | |
std::vector< std::array< double, 5 > > | generateIcosAxes (ProSHADE_comparePairwise *cmpObj, ProSHADE::ProSHADE_settings *settings, std::vector< std::array< double, 5 > > icosSymm, std::vector< std::array< double, 5 > > allCs, double axisErrorTolerance=0.1, int verbose=0) |
This function generates all the icosahedral symmetry elements from a pair detected by findIcosSymmetry function. More... | |
std::vector< std::array< double, 5 > > | generateIcosElements (std::vector< std::array< double, 5 > > symmAxes, ProSHADE::ProSHADE_settings *settings, int verbose=0) |
This function generates the 60 icosahedral symmetry group elements from the axes generated by generateIcosAxes function. More... | |
Public Attributes | |
std::vector< std::array< double, 5 > > | cnSymm |
This variable holds the complete Cyclic symmetry results. | |
std::vector< std::array< double, 5 > > | cnSymmClear |
This variable holds the gap corrected Cyclic symmetry results. | |
std::vector< std::vector< std::array< double, 6 > > > | dnSymm |
This variable holds the complete Dihedral symmetry results. | |
std::vector< std::vector< std::array< double, 6 > > > | dnSymmClear |
This variable holds the gap corrected Dihedral symmetry results. | |
std::vector< std::array< double, 5 > > | tetrAxes |
This variable holds the 7 unique axes of the tetrahedral symmetry, if they are found. | |
std::vector< std::array< double, 5 > > | octaAxes |
This variable holds the 13 unique axes of the octahedral symmetry, if they are found. | |
std::vector< std::array< double, 5 > > | icosAxes |
This variable holds the 31 unique axes of the octahedral symmetry, if they are found. | |
bool | inputStructureDataType |
This variable tells whether input data type is PDB or not. | |
Friends | |
class | ProSHADE_comparePairwise |
This class takes a single structure and the settings decided by the use and proceeds to compute all the requred information for symmetry detection in the constructor. Thus, once the object is constructed, the user only needs to call the accessor functions provided to get the deta and information rapidly.
Definition at line 48 of file ProSHADE_internal.h.
ProSHADE_internal::ProSHADE_symmetry::ProSHADE_symmetry | ( | ProSHADE::ProSHADE_settings * | settings | ) |
This is the constructor for the executive class which should be used by the user to obtain symmetry information about a single structure. This constructor contains all the symmetry running and reporting logic, while it relies on lower level functions for the specific computations. All computations are done immediately as the constructor is called, this may take some time, but all other operations are almost immediate and no more computations are required ever again. It should be used by the user instead of any other lower level classes.
[out] | X | Object containing all symmetry information that ProSHADE library can provide. |
Definition at line 1150 of file ProSHADE_internal.cpp.
ProSHADE_internal::ProSHADE_symmetry::ProSHADE_symmetry | ( | ) |
This is an empty constructor to allow access to the internal functions without the need for initialisation and computation of any results.
[out] | X | Empty object used only to allow access to internal functions without the need for initialisation. |
Definition at line 1133 of file ProSHADE_internal.cpp.
std::vector< std::array< double, 5 > > ProSHADE_internal::ProSHADE_symmetry::generateIcosAxes | ( | ProSHADE_comparePairwise * | cmpObj, |
ProSHADE::ProSHADE_settings * | settings, | ||
std::vector< std::array< double, 5 > > | icosSymms, | ||
std::vector< std::array< double, 5 > > | allCs, | ||
double | axisErrorTolerance = 0.1 , |
||
int | verbose = 0 |
||
) |
This function generates the 31 unique symmetry group elements for the icosahedral symmetry group. It does so by first determining the six C5 symmetries with cos(0.5) angle to each other, searching for missing ones if it cannot find them in the allCs list. It then searches for the ten C3 symmetries, which have the dihedral angle cos(-sqrt(5)/3) to the closer three C5s and the ( cos( 1.0 - ( sqrt(5)/3 ) ) ) angle to the further three C5s. If any are missing, it will attempt to locate these missing axes in the map. Finally, searching for the fifteen C2 symmetries with angles cos(0.0), cos(0.5) and cos(sqrt(3)/2) to the closest, mid-distance and furthest C5s. If all 31 axes are found, they are reported back, otherwise an empty vector is returned.
[in] | cmpObj | This is the object of ProSHADE_comparePairwise class, which allows access to the original inverse SO(3) FT map, so that missing axes detection can be done. |
[in] | settings | The ProSHADE_settings container class instance with information about how to write the HTML report. |
[in] | icosSymm | This is the pair of symmetry axes detected by the findIcosSymmetry function. |
[in] | allCs | The list of all symmetry axes in the whole object, from which the appropriate axes will be selected. |
[in] | axisErrorTolerance | The tolerance on the axis angle error. |
[in] | verbose | The amount of progress comments to be outputted to the std out. |
[out] | X | List of all elements making the icosahedral symmetry group, or empty if the function failed. |
Definition at line 4047 of file ProSHADE_SO3transform.cpp.
std::vector< std::array< double, 5 > > ProSHADE_internal::ProSHADE_symmetry::generateIcosElements | ( | std::vector< std::array< double, 5 > > | symmAxes, |
ProSHADE::ProSHADE_settings * | settings, | ||
int | verbose = 0 |
||
) |
[in] | symmAxes | The axes forming the tetrahedral symmetry as returned by the generateIcosAxes () function. |
[in] | settings | The ProSHADE_settings container class instance with information about how to write the HTML report. |
[in] | verbose | The amount of progress comments to be outputted to the std out. |
[out] | X | List of all elements making the icosahedral symmetry group. |
Definition at line 4626 of file ProSHADE_SO3transform.cpp.
std::vector< std::array< double, 5 > > ProSHADE_internal::ProSHADE_symmetry::generateOctaAxes | ( | ProSHADE_comparePairwise * | cmpObj, |
ProSHADE::ProSHADE_settings * | settings, | ||
std::vector< std::array< double, 5 > > | octaSymms, | ||
std::vector< std::array< double, 5 > > | allCs, | ||
double | axisErrorTolerance = 0.1 , |
||
int | verbose = 0 |
||
) |
This function searches the already detected Cn's for the presence of (cub)octahedral symmetry. It does so by first detecting three perpendicular C4 symmetries; in the case of one of these missing, it goes to search for it in the inverse SO(3) FT map. Then, it searches for the four C3 symmetries, which should have angle to all C4's of cos ( 1 / sqrt(3) ) - the dihedral angle. If not all C3 symmetries are found, it again goes to search for them. Finally, the function searches for the six C2 symmetries which should have angle of either cos(1/sqrt(2)) or cos(0) to the C4's. Again, if some are missing, it will search for them. The function then returns all 13 unique symmetry elements, or empty vector if it failed.
[in] | cmpObj | This is the object of ProSHADE_comparePairwise class, which allows access to the original inverse SO(3) FT map, so that missing axes detection can be done. |
[in] | settings | The ProSHADE_settings container class instance with information about how to write the HTML report. |
[in] | octaSymm | This is the pair of symmetry axes detected by the findOctaSymmetry function. |
[in] | allCs | The list of all symmetry axes in the whole object, from which the appropriate axes will be selected. |
[in] | axisErrorTolerance | The tolerance on the axis angle error. |
[in] | verbose | The amount of progress comments to be outputted to the std out. |
[out] | X | List of all elements making the (cub)octahedral symmetry group, or empty if the function failed. |
Definition at line 3339 of file ProSHADE_SO3transform.cpp.
std::vector< std::array< double, 5 > > ProSHADE_internal::ProSHADE_symmetry::generateOctaElements | ( | std::vector< std::array< double, 5 > > | symmAxes, |
ProSHADE::ProSHADE_settings * | settings, | ||
int | verbose = 0 |
||
) |
[in] | symmAxes | The axes forming the tetrahedral symmetry as returned by the generateOctaAxes () function. |
[in] | settings | The ProSHADE_settings container class instance with information about how to write the HTML report. |
[in] | verbose | The amount of progress comments to be outputted to the std out. |
[out] | X | List of all elements making the octahedral symmetry group. |
Definition at line 3921 of file ProSHADE_SO3transform.cpp.
std::vector< std::array< double, 5 > > ProSHADE_internal::ProSHADE_symmetry::generateTetrAxes | ( | ProSHADE_comparePairwise * | cmpObj, |
ProSHADE::ProSHADE_settings * | settings, | ||
std::vector< std::array< double, 5 > > | tetrSymms, | ||
std::vector< std::array< double, 5 > > | allCs, | ||
double | axisErrorTolerance = 0.1 , |
||
int | verbose = 0 |
||
) |
This function generates the tetrahedral symmetry group axes. It takes two C3 symmetries with angle acos(1/3) between them as a starting point and then proceeds to search for 2 most C3s, so that there would be alltogether four C3 symmetries, each with angle acos(1/3) to any other. If one of the C3's is not found, it will search for it in the inverse SO(3) FT map. Once these are located, three C2 symmetries perpendicular to each other are located each in between two C3s. Again, if there is one missing, the function will try to search for it. Once these are found, all the seven axes are returned; if not all seven axes are found, the function returns an empty vector.
[in] | cmpObj | This is the object of ProSHADE_comparePairwise class, which allows access to the original inverse SO(3) FT map, so that missing axes detection can be done. |
[in] | settings | The ProSHADE_settings container class instance with information about how to write the HTML report. |
[in] | tetrSymm | This is the pair of symmetry axes detected by the findTetrSymmetry function. |
[in] | allCs | The list of all symmetry axes in the whole object, from which the appropriate axes will be selected. |
[in] | axisErrorTolerance | The tolerance on the axis angle error. |
[in] | verbose | The amount of progress comments to be outputted to the std out. |
[out] | X | List of all axes making the tetrahedral symmetry group. |
Definition at line 2879 of file ProSHADE_SO3transform.cpp.
std::vector< std::array< double, 5 > > ProSHADE_internal::ProSHADE_symmetry::generateTetrElements | ( | std::vector< std::array< double, 5 > > | symmAxes, |
ProSHADE::ProSHADE_settings * | settings, | ||
int | verbose = 0 |
||
) |
[in] | symmAxes | The axes forming the tetrahedral symmetry as returned by the generateTetrAxes () function. |
[in] | settings | The ProSHADE_settings container class instance with information about how to write the HTML report. |
[in] | verbose | The amount of progress comments to be outputted to the std out. |
[out] | X | List of all elements making the tetrahedral symmetry group. |
Definition at line 3243 of file ProSHADE_SO3transform.cpp.
std::vector< std::array< double, 5 > > ProSHADE_internal::ProSHADE_symmetry::getCSymmetries | ( | ) |
This function gives a programmatical access to the results of the symmetry detection procedure C-symmetries list. This is a list of all detected C symmetries and is encoded as follows:
Each vector entry is a distinct C symmetry with its features encoded accordingly: [0] The fold of the symmetry [1] X-axis element of the axis for Angle-Axis representation [2] Y-axis element of the axis for Angle-Axis representation [3] Z-axis element of the axis for Angle-Axis representation [4] Average peak height for the symmetry
[out] | X | Vector of C-symmetries and their features. |
Definition at line 7241 of file ProSHADE_internal.cpp.
std::vector< std::vector< std::array< double, 6 > > > ProSHADE_internal::ProSHADE_symmetry::getDSymmetries | ( | ) |
This function gives a programmatical access to the results of the symmetry detection procedure D-symmetries list. This is a list of all detected D(ihedral) symmetries and is encoded as follows:
Each vector entry is a distinct D symmetry with its features encoded accordingly: [0] The fold of the symmetry [1] X-axis element average of the axis for Angle-Axis representation [2] Y-axis element average of the axis for Angle-Axis representation [3] Z-axis element average of the axis for Angle-Axis representation [4] <do not="" use>=""> [5] Average peak height for the symmetry
[out] | X | Vector of D-symmetries and their features. |
Definition at line 7266 of file ProSHADE_internal.cpp.
std::vector< std::array< double, 8 > > ProSHADE_internal::ProSHADE_symmetry::getRotFnPeaks | ( | ) |
This function gives a programmatical access to the results of the symmetry detection procedure peaks. Each peak reported in the output is where the overlay function of the object has value higher than the background as detected by the procedure. The output format is as follows:
Each vector entry is a distinct peak with its features encoded accordingly: [0] Euler Angle alpha value [1] Euler Angle beta value [2] Euler Angle gamma value [3] Peak height [4] X-axis element of the axis for Angle-Axis representation [5] Y-axis element of the axis for Angle-Axis representation [6] Z-axis element of the axis for Angle-Axis representation [7] Angle in radians for Angle-Axis representation
[out] | X | Vector of peaks and their features. |
Definition at line 7218 of file ProSHADE_internal.cpp.
void ProSHADE_internal::ProSHADE_symmetry::printResultsClear | ( | int | verbose | ) |
This version of the function prints the results to the screen with regard to the distribution of average peak height for the symmetries, so only the most relevant symmetries are shown. It has a chance of missing a true symmetry, but this typically happens in cases where the settings are disruptive for the prpcedure or when the computation has somehow failed.
[in] | verbose | Int value determining how verbal the function should be. It has currently no effect, though... |
Definition at line 5564 of file ProSHADE_internal.cpp.
void ProSHADE_internal::ProSHADE_symmetry::printResultsClearHTML | ( | ProSHADE::ProSHADE_settings * | settings | ) |
This version of the function prints the results to the screen with regard to the distribution of average peak height for the symmetries, so only the most relevant symmetries are shown. It has a chance of missing a true symmetry, but this typically happens in cases where the settings are disruptive for the procedure or when the computation has somehow failed.
[in] | settings | The ProSHADE_settings object containing all the information about the output and other things. |
Definition at line 5830 of file ProSHADE_internal.cpp.
void ProSHADE_internal::ProSHADE_symmetry::printResultsRequest | ( | std::string | symmetryType, |
unsigned int | symmetryFold, | ||
int | verbose | ||
) |
This version of the function prints the results to the screen with regard to the user specified expected symmetry. It not only searches for the requested symmetry in the already detected symmetries, but it also does extra searches if the requested symmetry cannot be found in the already detected ones.
[in] | symmetryType | String with the type of symmetry that is being sought after. Allowed values are C, D, T, O and I. |
[in] | symmetryFold | The fold of the symmetry being sought after. Only used for C and D symmetry types. |
[in] | verbose | Int value determining how verbal the function should be. It has currently no effect, though... |
Definition at line 3257 of file ProSHADE_internal.cpp.
void ProSHADE_internal::ProSHADE_symmetry::printResultsRequestHTML | ( | std::string | symmetryType, |
unsigned int | symmetryFold, | ||
ProSHADE::ProSHADE_settings * | settings | ||
) |
This version of the function prints the results to the HTML file with regard to the user specified expected symmetry. It not only searches for the requested symmetry in the already detected symmetries, but it also does extra searches if the requested symmetry cannot be found in the already detected ones.
[in] | symmetryType | String with the type of symmetry that is being sought after. Allowed values are C, D, T, O and I. |
[in] | symmetryFold | The fold of the symmetry being sought after. Only used for C and D symmetry types. |
[in] | settings | The ProSHADE_settings class instance containing all the information needed to write properly into the HTML output file. |
Definition at line 3661 of file ProSHADE_internal.cpp.