ProSHADE  0.7.5.1 (JAN 2021)
Protein Shape Detection
ProSHADE_distances.hpp
Go to the documentation of this file.
1 
23 //============================================ ProSHADE
25 
26 //============================================ Overinclusion protection
27 #ifndef __PROSHADE_DISTANCES__
28 #define __PROSHADE_DISTANCES__
29 
30 //============================================ ProSHADE_internal_distances Namespace
31 
40 {
43  bool isBandWithinShell ( proshade_unsign bandInQuestion, proshade_unsign shellInQuestion,
46  ProSHADE_settings* settings, proshade_unsign minCommonBands,
47  proshade_unsign minCommonShells, std::vector<proshade_double>* bandDists );
48  void allocateTrSigmaWorkspace ( proshade_unsign minSpheres, proshade_unsign intOrder, proshade_double*& obj1Vals, proshade_double*& obj2Vals,
49  proshade_double*& GLabscissas, proshade_double*& glWeights, proshade_complex*& radiiVals );
50  void computeSphericalHarmonicsMagnitude ( ProSHADE_internal_data::ProSHADE_data* obj, proshade_unsign band, proshade_unsign order,
51  proshade_unsign radius, proshade_double* result );
53  proshade_unsign bandIter, proshade_unsign orderIter, proshade_complex* radiiVals, proshade_unsign integOrder,
54  proshade_double* abscissas, proshade_double* weights, proshade_double integRange, proshade_double sphereDist );
56  proshade_unsign bandIter, proshade_unsign orderIter, proshade_double* obj1Vals, proshade_double* obj2Vals,
57  proshade_unsign integOrder, proshade_double* abscissas, proshade_double* weights, proshade_double sphereDist );
58  void releaseTrSigmaWorkspace ( proshade_double*& obj1Vals, proshade_double*& obj2Vals, proshade_double*& GLabscissas,
59  proshade_double*& glWeights, proshade_complex*& radiiVals );
61  ProSHADE_settings* settings );
63  ProSHADE_settings* settings );
65  ProSHADE_settings* settings );
67  ProSHADE_settings* settings );
68  void allocateInvSOFTWorkspaces ( proshade_complex*& work1, proshade_complex*& work2, proshade_double*& work3, proshade_unsign band );
69  void prepareInvSOFTPlan ( fftw_plan* inverseSO3, proshade_unsign band, fftw_complex* work1, proshade_complex* invCoeffs );
70  void releaseInvSOFTMemory ( proshade_complex*& work1, proshade_complex*& work2, proshade_double*& work3 );
72  ProSHADE_settings* settings );
74  ProSHADE_settings* settings );
75 }
76 
77 #endif
ProSHADE_internal_distances::normaliseEMatrices
void normaliseEMatrices(ProSHADE_internal_data::ProSHADE_data *obj1, ProSHADE_internal_data::ProSHADE_data *obj2, ProSHADE_settings *settings)
This function normalises the E matrices.
Definition: ProSHADE_distances.cpp:577
ProSHADE_internal_distances::allocateTrSigmaWorkspace
void allocateTrSigmaWorkspace(proshade_unsign minSpheres, proshade_unsign intOrder, proshade_double *&obj1Vals, proshade_double *&obj2Vals, proshade_double *&GLabscissas, proshade_double *&glWeights, proshade_complex *&radiiVals)
This helper function is responsible for allocating the workspace memory required for trace sigma desc...
Definition: ProSHADE_distances.cpp:317
ProSHADE_internal_distances::computeWeightsForEMatricesForLM
proshade_double computeWeightsForEMatricesForLM(ProSHADE_internal_data::ProSHADE_data *obj1, ProSHADE_internal_data::ProSHADE_data *obj2, proshade_unsign bandIter, proshade_unsign orderIter, proshade_double *obj1Vals, proshade_double *obj2Vals, proshade_unsign integOrder, proshade_double *abscissas, proshade_double *weights, proshade_double sphereDist)
This function computes the E matrix weight values for a given band and order and saves these into the...
Definition: ProSHADE_distances.cpp:438
ProSHADE_internal_data::ProSHADE_data
This class contains all inputed and derived data for a single structure.
Definition: ProSHADE_data.hpp:49
ProSHADE_internal_distances::releaseInvSOFTMemory
void releaseInvSOFTMemory(proshade_complex *&work1, proshade_complex *&work2, proshade_double *&work3)
This function releases the memory used for computation of the inverse SOFT transform.
Definition: ProSHADE_distances.cpp:837
ProSHADE_internal_distances::isBandWithinShell
bool isBandWithinShell(proshade_unsign bandInQuestion, proshade_unsign shellInQuestion, ProSHADE_internal_spheres::ProSHADE_sphere **spheres)
This function checks if a band is available for a given shell.
Definition: ProSHADE_distances.cpp:143
ProSHADE_internal_spheres::ProSHADE_sphere
This class contains all inputed and derived data for a single sphere.
Definition: ProSHADE_spheres.hpp:49
ProSHADE_internal_distances::computeRotationunctionDescriptor
proshade_double computeRotationunctionDescriptor(ProSHADE_internal_data::ProSHADE_data *obj1, ProSHADE_internal_data::ProSHADE_data *obj2, ProSHADE_settings *settings)
This function computes the rotation function descriptor value between two objects.
Definition: ProSHADE_distances.cpp:911
ProSHADE_internal_distances
This namespace contains the functions used for computing distances between two structures as represen...
ProSHADE_internal_distances::computeRRPPearsonCoefficients
void computeRRPPearsonCoefficients(ProSHADE_internal_data::ProSHADE_data *obj1, ProSHADE_internal_data::ProSHADE_data *obj2, ProSHADE_settings *settings, proshade_unsign minCommonBands, proshade_unsign minCommonShells, std::vector< proshade_double > *bandDists)
This function gets the Pearson's coefficients or all bands between two objects.
Definition: ProSHADE_distances.cpp:215
ProSHADE_internal_distances::prepareInvSOFTPlan
void prepareInvSOFTPlan(fftw_plan *inverseSO3, proshade_unsign band, fftw_complex *work1, proshade_complex *invCoeffs)
This function prepares the FFTW plan for the inverse SO(3) transform.
Definition: ProSHADE_distances.cpp:790
ProSHADE_settings
This class stores all the settings and is passed to the executive classes instead of a multitude of p...
Definition: ProSHADE_settings.hpp:86
ProSHADE_internal_distances::computeEMatrices
void computeEMatrices(ProSHADE_internal_data::ProSHADE_data *obj1, ProSHADE_internal_data::ProSHADE_data *obj2, ProSHADE_settings *settings)
This function computes the complete E matrices and their weights between any two objects.
Definition: ProSHADE_distances.cpp:515
ProSHADE_internal_distances::computeInverseSOFTTransform
void computeInverseSOFTTransform(ProSHADE_internal_data::ProSHADE_data *obj1, ProSHADE_internal_data::ProSHADE_data *obj2, ProSHADE_settings *settings)
This function computes the inverse SO(3) transform.
Definition: ProSHADE_distances.cpp:859
ProSHADE_internal_distances::computeEMatricesForLM
void computeEMatricesForLM(ProSHADE_internal_data::ProSHADE_data *obj1, ProSHADE_internal_data::ProSHADE_data *obj2, proshade_unsign bandIter, proshade_unsign orderIter, proshade_complex *radiiVals, proshade_unsign integOrder, proshade_double *abscissas, proshade_double *weights, proshade_double integRange, proshade_double sphereDist)
This function computes the E matrix un-weighted values for a given band and order and saves these int...
Definition: ProSHADE_distances.cpp:375
ProSHADE_internal_distances::computeSphericalHarmonicsMagnitude
void computeSphericalHarmonicsMagnitude(ProSHADE_internal_data::ProSHADE_data *obj, proshade_unsign band, proshade_unsign order, proshade_unsign radius, proshade_double *result)
This function computes the magnitude of a particular spherical harmonics position for a given object,...
Definition: ProSHADE_distances.cpp:346
ProSHADE_internal_distances::computeTraceSigmaDescriptor
proshade_double computeTraceSigmaDescriptor(ProSHADE_internal_data::ProSHADE_data *obj1, ProSHADE_internal_data::ProSHADE_data *obj2, ProSHADE_settings *settings)
This function computes the trace sigma descriptor value between two objects.
Definition: ProSHADE_distances.cpp:621
ProSHADE_internal_distances::releaseTrSigmaWorkspace
void releaseTrSigmaWorkspace(proshade_double *&obj1Vals, proshade_double *&obj2Vals, proshade_double *&GLabscissas, proshade_double *&glWeights, proshade_complex *&radiiVals)
This helper function is responsible for deleting the workspace memory required for trace sigma descri...
Definition: ProSHADE_distances.cpp:483
ProSHADE_internal_distances::allocateInvSOFTWorkspaces
void allocateInvSOFTWorkspaces(proshade_complex *&work1, proshade_complex *&work2, proshade_double *&work3, proshade_unsign band)
This function allocates the workspaces required to compute the inverse SOFT transform.
Definition: ProSHADE_distances.cpp:766
ProSHADE_internal_distances::computeEnergyLevelsDescriptor
proshade_double computeEnergyLevelsDescriptor(ProSHADE_internal_data::ProSHADE_data *obj1, ProSHADE_internal_data::ProSHADE_data *obj2, ProSHADE_settings *settings)
This function computes the energy levels descriptor value between two objects.
Definition: ProSHADE_distances.cpp:165
ProSHADE_wignerMatrices.hpp
This header declares the functions required to compute the Wigner D matrices.
ProSHADE_internal_distances::generateSO3CoeffsFromEMatrices
void generateSO3CoeffsFromEMatrices(ProSHADE_internal_data::ProSHADE_data *obj1, ProSHADE_internal_data::ProSHADE_data *obj2, ProSHADE_settings *settings)
This function converts the E matrices to SO(3) coefficients.
Definition: ProSHADE_distances.cpp:706