ProSHADE  0.6.5 (NOV 2018)
Protein Shape Descriptors and Symmetry Detection
ProSHADE_internal::ProSHADE_comparePairwise Class Reference

This is the executive class responsible for comparing strictly two structures. More...

#include <ProSHADE_internal.h>

Public Member Functions

 ProSHADE_comparePairwise (ProSHADE_data *cmpObj1, ProSHADE_data *cmpObj2, double mPower, std::vector< int > ignoreL, unsigned int order, ProSHADE::ProSHADE_settings *settings)
 Contructor for the ProSHADE_comparePairwise class. More...
 
 ~ProSHADE_comparePairwise ()
 Destructor for the ProSHADE_comparePairwise class. More...
 
void precomputeTrSigmaDescriptor ()
 This function computes the E matrices required for the trace sigma descriptor, the rotation function descriptor and SO3 transform. More...
 
void getSO3InverseMap (ProSHADE::ProSHADE_settings *settings)
 This function is responsible for computing the SO3 inverse transform. More...
 
std::array< double, 3 > getEulerAngles (ProSHADE::ProSHADE_settings *settings, double *correlation=nullptr)
 This function finds the highest peak in the SO3 inverse transform map and sets it as the optimal overlay angle. More...
 
void setEulerAngles (double alpha, double beta, double gamma)
 This function is used to set the Euler angles for further processing. More...
 
std::vector< std::array< double, 8 > > getSO3Peaks (ProSHADE::ProSHADE_settings *settings, double noIQRs=1.5, bool freeMem=true, int peakSize=1, double realDist=0.4, int verbose=1)
 This function detects 'true' peaks from the background of the SO3 inverse transform map. More...
 
bool checkPeakExistence (double alpha, double beta, double gamma, ProSHADE::ProSHADE_settings *settings, double passThreshold=0.3)
 This function checks if a peak Euler angles actually produce something the descriptors would agree is similar. More...
 
std::vector< std::array< double, 5 > > findCnSymmetry (std::vector< std::array< double, 8 > > peaks, ProSHADE::ProSHADE_settings *settings, double axisErrorTolerance=0.0, bool freeMem=true, double percAllowedToMiss=0.33, int verbose=1)
 This function attempts to detect the C-symmetries present in the list of overlay peaks. More...
 
std::vector< std::array< double, 5 > > findCnSymmetryClear (std::vector< std::array< double, 5 > > CnSymm, ProSHADE::ProSHADE_settings *settings, double maxGap=0.2, bool *pf=nullptr)
 This function takes the detected C-symmetries list, removes low probability symmetries and sorts as to make the output clearer. More...
 
std::vector< std::vector< std::array< double, 6 > > > findDnSymmetry (std::vector< std::array< double, 5 > > CnSymm, double axisErrorTolerance=0.1)
 This function detects dihedral (D) symmetries from the list of already detected C symmetries. More...
 
std::vector< std::vector< std::array< double, 6 > > > findDnSymmetryClear (std::vector< std::vector< std::array< double, 6 > > > DnSymm, ProSHADE::ProSHADE_settings *settings, double maxGap=0.2, bool *pf=nullptr)
 This function sorts the D symmetries and removed low probability ones. More...
 
std::vector< std::vector< std::array< double, 5 > > > findTetrSymmetry (std::vector< std::array< double, 5 > > CnSymm, double *tetrSymmPeakAvg, double axisErrorTolerance=0.1)
 This function detects the tetrahedral symmetries in the structure. More...
 
std::vector< std::vector< std::array< double, 5 > > > findOctaSymmetry (std::vector< std::array< double, 5 > > CnSymm, double *octaSymmPeakAvg, double axisErrorTolerance=0.1)
 This function detects the octahedral symmetries in the structure. More...
 
std::vector< std::vector< std::array< double, 5 > > > findIcosSymmetry (std::vector< std::array< double, 5 > > CnSymm, double *icosSymmPeakAvg, double axisErrorTolerance=0.1)
 This function detects the icosahedral symmetries in the structure. More...
 
void generateWignerMatrices (ProSHADE::ProSHADE_settings *settings)
 This function is responsible for computing the Wigner D matrices for symmetry detection purposes. More...
 
double getRotCoeffDistance (ProSHADE::ProSHADE_settings *settings)
 This function computes the full rotation function descriptor distances. More...
 
double maxAvgPeakForSymmetry (double X, double Y, double Z, double Angle, ProSHADE::ProSHADE_settings *settings)
 This function takes angle and axis and checks the SO(3) map for this specific symmetry only. More...
 
void freeInvMap ()
 This function frees the SOFT inverse map. More...
 
double getPeakThreshold ()
 This is an accessor function to the peak height threshold. More...
 
void rotateStructure (ProSHADE_data *cmpObj1, ProSHADE::ProSHADE_settings *settings, std::string saveName, int verbose=0, std::string axOrd="xyz", bool internalUse=false)
 Rotates the density map be given Angle-Axis rotation using the Wigner matrices and inverse spharical harmonics transform. More...
 
std::array< double, 4 > getTranslationFunctionMap (ProSHADE_data *obj1, ProSHADE_data *obj2, double *ob2XMov=nullptr, double *ob2YMov=nullptr, double *ob2ZMov=nullptr)
 Computes the translation function for the two objects and returns the position of the highest peak. More...
 

Detailed Description

This class is used to compare strictly two structures and is currently used for symmetry detection, but should not really be used for the distance computation, as there is designated class for this (ProSHADE_distances) and this class uses the current and best tested approach.

Definition at line 358 of file ProSHADE_internal.h.

Constructor & Destructor Documentation

§ ProSHADE_comparePairwise()

ProSHADE_internal::ProSHADE_comparePairwise::ProSHADE_comparePairwise ( ProSHADE_data cmpObj1,
ProSHADE_data cmpObj2,
double  mPower,
std::vector< int >  ignoreL,
unsigned int  order,
ProSHADE::ProSHADE_settings settings 
)

This is the constructor for the executive class responsible for comparing two structures with identical settings. This is currently used only for symmetry detection as the distance computation has its designated class ProSHADE_distances. The constructor is responsible for sanity checks, setting up parameters and determining the integration variables.

Parameters
[in]cmpObj1This is the pointer to the first structure data holding object.
[in]cmpObj2This is the pointer to the second structure data holding object.
[in]mPowerThis is the weighting parameter for energy level distances computation.
[in]ignoreLThis vector should contain all bands to be ignored in the computations.
[in]orderThis is the order of the Gauss-Legendre integration to be used.
[in]settingsProSHADE_settings container class with information about the HTML report printing, if needed.
[out]XObject capable of comparing the two structures.
Warning
This class should not be directly accessed by the user.

Definition at line 502 of file ProSHADE_internal.cpp.

§ ~ProSHADE_comparePairwise()

ProSHADE_internal::ProSHADE_comparePairwise::~ProSHADE_comparePairwise ( )

The simple destructor is responsible for checking that all the memory has been appropriately released before and if not, then it releases it itself.

Warning
This is an internal function which should not be used by the user.

Definition at line 604 of file ProSHADE_internal.cpp.

Member Function Documentation

§ checkPeakExistence()

bool ProSHADE_internal::ProSHADE_comparePairwise::checkPeakExistence ( double  alpha,
double  beta,
double  gamma,
ProSHADE::ProSHADE_settings settings,
double  passThreshold = 0.3 
)

This function is used internally to check the passing 'true' peaks for actually leading to structure overlay which is similar to the original. This is done in attempt to remove 'false' peaks where the background noise produces a peak which however is not symmetry related.

Parameters
[in]alphaThe Euler ZXZ alpha angle of the peak to be checked.
[in]betaThe Euler ZXZ beta angle of the peak to be checked.
[in]gammaThe Euler ZXZ gamma angle of the peak to be checked.
[in]settingsThe Settings class with information about the HTML report writing, if needed.
[in]passThresholdThe distance threshold for passing as 'true' peak instead of being considered 'background' peak.
[out]XBoolean value whether the peak is 'true'.
Warning
This is an internal function which should not be used by the user.

Definition at line 2594 of file ProSHADE_SO3transform.cpp.

§ findCnSymmetry()

std::vector< std::array< double, 5 > > ProSHADE_internal::ProSHADE_comparePairwise::findCnSymmetry ( std::vector< std::array< double, 8 > >  peaks,
ProSHADE::ProSHADE_settings settings,
double  axisErrorTolerance = 0.0,
bool  freeMem = true,
double  percAllowedToMiss = 0.33,
int  verbose = 1 
)

This function takes the list of SO3 inverse transform peaks (courtesy of getSO3Peaks) and searches for C symmetries. This is done by first grouping the peaks by the axis of rotation and then for each group with the same axis, repetitively finding the smallest angle between the angle values and checking if such symmetry is supported by the rest of the peaks.

Parameters
[in]peaksThe list of SO3 inverse transform peaks as produced by the getSO3Peaks function.
[in]settingsThe ProSHADE_settings container class instance with instructions on how to write the HTML report.
[in]axisErrorToleranceTolerance value on matching the peaks axis together.
[in]freeMemShould the memory be released, or is this search going to be repeated?
[in]percAllowedToMissThe percentage (0 to 1) of peaks that are allowed to be missing.
[in]verboseAn integer value stating how quiet (0) or loud (4) the function reports should be. In this case, it will be quiet completely unless the value is 4.
[out]XList of C symmetries detected.
Warning
This is an internal function which should not be used by the user.

Definition at line 1292 of file ProSHADE_SO3transform.cpp.

§ findCnSymmetryClear()

std::vector< std::array< double, 5 > > ProSHADE_internal::ProSHADE_comparePairwise::findCnSymmetryClear ( std::vector< std::array< double, 5 > >  CnSymm,
ProSHADE::ProSHADE_settings settings,
double  maxGap = 0.2,
bool *  pf = nullptr 
)

This is more of a visualisation/presentation function, as it does not compute anything important, it just removes C symmetries which have considerably higher one already present and sorts as to make the highest symmetry order first (if several symmetries are found) or according to the peak height (if there are many reported symmetries and the user will have to interfere).

Parameters
[in]CnSymmList of the detected C symmetries as produced by the findCnSymmetry function.
[in]settingsThe ProSHADE_settings container class instance with information about how to write the HTML report.
[in]maxGapThe minimum gap required between the symmetries to cut the output.
[in]pfPointer to boolean variable to decide whether full or clear output should be printed.
[out]XList of C symmetries passing these more stringent conditions.
Warning
This is an internal function which should not be used by the user.

Definition at line 2128 of file ProSHADE_SO3transform.cpp.

§ findDnSymmetry()

std::vector< std::vector< std::array< double, 6 > > > ProSHADE_internal::ProSHADE_comparePairwise::findDnSymmetry ( std::vector< std::array< double, 5 > >  CnSymm,
double  axisErrorTolerance = 0.1 
)

Given that the C symmetries have been determined, it is now rather simple to detect the D symmetries, as there is a D symmetry whenever two C symmetries have dot product of their axes 0 (are perpendicular). Thus, this function does all the required work and reports the detected D-symmetries.

Parameters
[in]CnSymmList of the detected C symmetries as produced by the findCnSymmetry function.
[in]axisErrorToleranceTolerance value on matching the peaks axis together.
[out]XList of detected D symmetries.
Warning
This is an internal function which should not be used by the user.

Definition at line 2211 of file ProSHADE_SO3transform.cpp.

§ findDnSymmetryClear()

std::vector< std::vector< std::array< double, 6 > > > ProSHADE_internal::ProSHADE_comparePairwise::findDnSymmetryClear ( std::vector< std::vector< std::array< double, 6 > > >  DnSymm,
ProSHADE::ProSHADE_settings settings,
double  maxGap = 0.2,
bool *  pf = nullptr 
)

Given that the D symmetries are determined, they are now sorted by the symmetry order and peak height. Also, D symmetries with are much lower that the top ones are omitted to make for a clearer results screen. Finally, only the highest two D symmetries are reported and the Dx-x instead of all matching.

Parameters
[in]CnSymmList of detected D symmetries as produced by the findDnSymmetry function.
[in]settingsThe ProSHADE_settings container class instance with information about how to write the HTML report.
[in]maxGapThe minimum gap required between the symmetries to cut the output.
[in]pfPointer to boolean variable to decide whether full or clear output should be printed.
[out]XList of D symmetries passing these more stringent conditions.
Warning
This is an internal function which should not be used by the user.

Definition at line 2418 of file ProSHADE_SO3transform.cpp.

§ findIcosSymmetry()

std::vector< std::vector< std::array< double, 5 > > > ProSHADE_internal::ProSHADE_comparePairwise::findIcosSymmetry ( std::vector< std::array< double, 5 > >  CnSymm,
double *  icosSymmPeakAvg,
double  axisErrorTolerance = 0.1 
)

This function takes the C symmetries list and checks for existence of C5 symmetries. All of these are then tested against all the C3 symmetries found in the data and if the relative angle of their axes is equal to acos ( sqrt ( 5.0 ) / 3.0 ), then the structure has to have the icosahedral symmetry.

Parameters
[in]CnSymmList of the already detected C symmetries.
[in]icosSymmPeakAvgA pointer to variable holding the average peak value for all icosahedral symmetry supporting C symmetries.
[in]axisErrorToleranceThe tolerance on the axes angle error.
[out]XList of icosahedral symmetry showing C symmetries.
Warning
This is an internal function which should not be used by the user.

Definition at line 2795 of file ProSHADE_SO3transform.cpp.

§ findOctaSymmetry()

std::vector< std::vector< std::array< double, 5 > > > ProSHADE_internal::ProSHADE_comparePairwise::findOctaSymmetry ( std::vector< std::array< double, 5 > >  CnSymm,
double *  octaSymmPeakAvg,
double  axisErrorTolerance = 0.1 
)

This function takes the C symmetries list and checks for existence of C4 symmetries. All of these are then tested against all the different C3 symmetries found in the data and if the relative angle of their axes is equal to PI - acos ( 1.0 / 3.0 ), then the structure has to have the octahedral symmetry.

Parameters
[in]CnSymmList of the already detected C symmetries.
[in]octaSymmPeakAvgA pointer to variable holding the average peak value for all octahedral symmetry supporting C symmetries.
[in]axisErrorToleranceThe tolerance on the axes angle error.
[out]XList of octahedral symmetry showing C symmetries.
Warning
This is an internal function which should not be used by the user.

Definition at line 2715 of file ProSHADE_SO3transform.cpp.

§ findTetrSymmetry()

std::vector< std::vector< std::array< double, 5 > > > ProSHADE_internal::ProSHADE_comparePairwise::findTetrSymmetry ( std::vector< std::array< double, 5 > >  CnSymm,
double *  tetrSymmPeakAvg,
double  axisErrorTolerance = 0.1 
)

This function takes the C symmetries list and checks for existence of C3 symmetries. All of these are then tested against all the different C3 symmetries found in the data and if the relative angle of their axes is equal to acos ( 1.0 / 3.0 ), then the structure has to have the tetrahedral symmetry.

Parameters
[in]CnSymmList of the already detected C symmetries.
[in]tetrSymmPeakAvgA pointer to variable holding the average peak value for all tetrahedral symmetry supporting C symmetries.
[in]axisErrorToleranceThe tolerance on the axes angle error.
[out]XList of tetrahedral symmetry showing C symmetries.
Warning
This is an internal function which should not be used by the user.

Definition at line 2633 of file ProSHADE_SO3transform.cpp.

§ freeInvMap()

void ProSHADE_internal::ProSHADE_comparePairwise::freeInvMap ( )

The need for this function is specific to the symmetry detection, where the inverse SO(3) Fourier transform map is required to linger a bit longer than just the peak detection. This is because missing axis searches do require it, so that specific timing of the memory freeing needs to be done.

Warning
This is an internal function which should not be used by the user.

Definition at line 2084 of file ProSHADE_SO3transform.cpp.

§ generateWignerMatrices()

void ProSHADE_internal::ProSHADE_comparePairwise::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 symmetry detection class. The main purpose of the Wigner D matrices is to rotate the spherical harmonics coefficients.

Parameters
[in]settingsThe ProSHADE_settings container class with information about how to write the HTML report.
Warning
This is an internal function which should not be used by the user.

Definition at line 58 of file ProSHADE_wigner.cpp.

§ getEulerAngles()

std::array< double, 3 > ProSHADE_internal::ProSHADE_comparePairwise::getEulerAngles ( ProSHADE::ProSHADE_settings settings,
double *  correlation = nullptr 
)

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.

Parameters
[in]settingsThe settings class container with information about how to print the HTML report.
[in]correlationThe maximum peak height for the detected best overlay angles.
[out]XThis 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.
Warning
This is an internal function which should not be used by the user.

Definition at line 429 of file ProSHADE_SO3transform.cpp.

§ getPeakThreshold()

double ProSHADE_internal::ProSHADE_comparePairwise::getPeakThreshold ( )

This function provides access to the private class member peak threshold. It is needed for symmetry detection, where missing axes needs to know the minimum height required, but does not have access to the private members.

Parameters
[out]XThe value of the peak threshold used in peak detection.
Warning
This is an internal function which should not be used by the user.

Definition at line 2108 of file ProSHADE_SO3transform.cpp.

§ getRotCoeffDistance()

double ProSHADE_internal::ProSHADE_comparePairwise::getRotCoeffDistance ( 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.

Parameters
[in]settingsProSHADE_settings container class with information about the HTML report printing, if needed.
[out]XA double precision number which is the distance.
Warning
This is an internal function which should not be used by the user.

Definition at line 988 of file ProSHADE_descriptors.cpp.

§ getSO3InverseMap()

void ProSHADE_internal::ProSHADE_comparePairwise::getSO3InverseMap ( ProSHADE::ProSHADE_settings settings)

The SO3 transform and its inverse are required by the symmetry detection class to obtain peaks where the structure overlays itself well and this function provides the transform to allow it.

Parameters
[in]settingsThe settings class container with information about how to print the HTML report.
Warning
This is an internal function which should not be used by the user.

Definition at line 58 of file ProSHADE_SO3transform.cpp.

§ getSO3Peaks()

std::vector< std::array< double, 8 > > ProSHADE_internal::ProSHADE_comparePairwise::getSO3Peaks ( ProSHADE::ProSHADE_settings settings,
double  noIQRs = 1.5,
bool  freeMem = true,
int  peakSize = 1,
double  realDist = 0.4,
int  verbose = 1 
)

This function looks through the SO3 inverse map and searches for peaks. This is done by first finding all peaks (i.e. values which are higher than all neighouring values) and then determining background threshold so that only 'true' peaks are processed. It then optimises the 'true' peaks and converts them onto the angle-axis representation. This if finally returned by the function.

Parameters
[in]settingsThe ProSHADE_settings class container with infromation about the HTML report output, if needed.
[in]noIQRsThe number of interquartile ranges to determine the peak threshold.
[in]freeMemShould the memory be released, or will there be more peak operations to follow?
[in]peakSizeHow many neighbours in each direction need to be less than the cental value to be considered a peak?
[in]realDistThe full rotation function descriptor value which needs to be surpassed by a 'missing' peak in order to be 'found'.
[in]verboseInt value specifying how quiet (0) or loud (4) the function should be when reporting progress.
[out]XA list of peaks with high overlay value.
Warning
This is an internal function which should not be used by the user.

Definition at line 621 of file ProSHADE_SO3transform.cpp.

§ getTranslationFunctionMap()

std::array< double, 4 > ProSHADE_internal::ProSHADE_comparePairwise::getTranslationFunctionMap ( ProSHADE_data obj1,
ProSHADE_data obj2,
double *  ob2XMov = nullptr,
double *  ob2YMov = nullptr,
double *  ob2ZMov = nullptr 
)

This function takes two ProSHADE_data objects and proceeds to make sure they have the same dimensions and map sampling. This does change the maps, so keep this in mind. Then, it computes Fourier transform on both structures, combines the Fourier coefficients which are now on the same scale and inverts the Fourier transform. This results in the translation function map. Finally, finding the highest peak and its coordinates, the optimal overlay translation is obtained and transformed into translation in angstroms instead of indices. This translation in ANGSTROMS is then returned.

Parameters
[in]obj1Pointer to the ProSHADE_data object which has the density map to be matched against obj2.
[in]obj2Pointer to the ProSHADE_data object which has the density map to be matched against obj1.
[out]XArray of three doubles specifying the optimal overlay translation in angstroms between the two input objects.
Warning
This is an internal function which should not be used by the user.

Definition at line 4754 of file ProSHADE_SO3transform.cpp.

§ maxAvgPeakForSymmetry()

double ProSHADE_internal::ProSHADE_comparePairwise::maxAvgPeakForSymmetry ( double  X,
double  Y,
double  Z,
double  Angle,
ProSHADE::ProSHADE_settings settings 
)

This function takes angle and axis and checks the SO(3) map for this specific symmetry only.

Parameters
[in]XThe X-axis element of the angle-axis rotation representation for which to check.
[in]YThe Y-axis element of the angle-axis rotation representation for which to check.
[in]ZThe Z-axis element of the angle-axis rotation representation for which to check.
[in]AngleThe angle of the angle-axis rotation representation for which to check.
[in]settingsThe ProSHADE_settings container class instance with instructions on how to write the HTML report.
[out]XValue of the average peak height for the best possible axis found.
Warning
This is an internal function which should not be used by the user.

Definition at line 1143 of file ProSHADE_SO3transform.cpp.

§ precomputeTrSigmaDescriptor()

void ProSHADE_internal::ProSHADE_comparePairwise::precomputeTrSigmaDescriptor ( )

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 symmetry detection.

Warning
This is an internal function which should not be used by the user.

Definition at line 454 of file ProSHADE_descriptors.cpp.

§ rotateStructure()

void ProSHADE_internal::ProSHADE_comparePairwise::rotateStructure ( ProSHADE_data cmpObj1,
ProSHADE::ProSHADE_settings settings,
std::string  saveName,
int  verbose = 0,
std::string  axOrd = "xyz",
bool  internalUse = false 
)

This function takes the density map of a structure, computes the Wigner D matrices using the Euler angles already set for the structure (this needs to be done beforehand) and then proceeds to multiply the spherical coordiantes by the Wigner D matrices, thus computing the rotation in the spherical harmonics (Fourier) space. Consequently, it inverts the rotated spherical harmonics coefficients and coverts the resulting sphere mapped values back to Cartesian grid, which is finally outputted as a map. The result is a density map rotated by the required Euler angles without any interpolation in the real space, but with interpolation to and from spheres.

Parameters
[in]cmpObj1Pointer to the ProSHADE_data object which has the density map and angles set.
[in]settingsThe ProSHADE_settings container class instance with information about how to write the HTML report.
[in]saveNameString of the name to where the rotated structure should ba outputted to.
[in]verboseInt specifying how loud the function should be when reporting progress.
[in]axOrdString of three letters, specifying the order of axes, allowed characters are 'x', 'y' and 'z' in any order.
[in]internalUseIf this is for internal use, save some time on 0, 0, 0 angles.
[in]noMapBoolean value specifying whether a map is at all needed and possible to obtain.
Warning
This is an internal function which should not be used by the user.

Definition at line 5641 of file ProSHADE_SO3transform.cpp.

§ setEulerAngles()

void ProSHADE_internal::ProSHADE_comparePairwise::setEulerAngles ( double  alpha,
double  beta,
double  gamma 
)

This function is called when the proper rotation between two structures has already been determined to set this rotation for further processing, specifically computation of Wigner D matrices using these. Except for setting the internal variables, it makes sure the sanity checks further down the processing pipeline are passed, so it needs to be used instead of setting the internal variables directly.

Parameters
[in]alphaThe alpha angle of the optimal overlay Euler angles set.
[in]betaThe beta angle of the optimal overlay Euler angles set.
[in]gammaThe gamma angle of the optimal overlay Euler angles set.
Warning
This is an internal function which should not be used by the user.

Definition at line 401 of file ProSHADE_SO3transform.cpp.


The documentation for this class was generated from the following files: