ProSHADE  0.7.5.4 (MAR 2021)
Protein Shape Detection
ProSHADE_overlay.hpp
Go to the documentation of this file.
1 
22 //==================================================== ProSHADE
23 #include "ProSHADE_symmetry.hpp"
24 
25 //==================================================== Overinclusion protection
26 #ifndef __PROSHADE_OVERLAY__
27 #define __PROSHADE_OVERLAY__
28 
29 //==================================================== ProSHADE_internal_overlay Namespace
37 {
39  ProSHADE_internal_data::ProSHADE_data* movingStructure, proshade_double* eulA, proshade_double* eulB,
40  proshade_double* eulG );
42  ProSHADE_internal_data::ProSHADE_data* movingStructure, proshade_double* trsX,
43  proshade_double* trsY, proshade_double* trsZ, proshade_double eulA, proshade_double eulB,
44  proshade_double eulG );
45  void computeBeforeAfterZeroCounts ( proshade_unsign* addXPre, proshade_unsign* addYPre, proshade_unsign* addZPre, proshade_unsign* addXPost,
46  proshade_unsign* addYPost, proshade_unsign* addZPost, proshade_unsign xDim, proshade_unsign yDim,
47  proshade_unsign zDim, proshade_unsign xDimIndices, proshade_unsign yDimIndices, proshade_unsign zDimIndices );
48  void paddMapWithZeroes ( proshade_double* oldMap, proshade_double*& newMap, proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim,
49  proshade_unsign xDimIndices, proshade_unsign yDimIndices, proshade_unsign zDimIndices, proshade_unsign addXPre,
50  proshade_unsign addYPre, proshade_unsign addZPre );
51  void allocateTranslationFunctionMemory ( fftw_complex*& tmpIn1, fftw_complex*& tmpOut1, fftw_complex*& tmpIn2, fftw_complex*& tmpOut2, fftw_complex*& resIn,
52  fftw_complex*& resOut, fftw_plan& forwardFourierObj1, fftw_plan& forwardFourierObj2, fftw_plan& inverseFourierCombo,
53  proshade_unsign xD, proshade_unsign yD, proshade_unsign zD );
54  void combineFourierForTranslation ( fftw_complex* tmpOut1, fftw_complex* tmpOut2, fftw_complex*& resOut, proshade_unsign xD, proshade_unsign yD,
55  proshade_unsign zD );
56  void findHighestValueInMap ( fftw_complex* resIn, proshade_unsign xD, proshade_unsign yD, proshade_unsign zD, proshade_double* trsX,
57  proshade_double* trsY, proshade_double* trsZ, proshade_double* mapPeak );
58  void freeTranslationFunctionMemory ( fftw_complex*& tmpIn1, fftw_complex*& tmpOut1, fftw_complex*& tmpIn2, fftw_complex*& tmpOut2, fftw_complex*& resOut,
59  fftw_plan& forwardFourierObj1, fftw_plan& forwardFourierObj2, fftw_plan& inverseFourierCombo );
60  void computeAngularThreshold ( std::vector<proshade_double>* lonCO, std::vector<proshade_double>* latCO, proshade_unsign angRes );
61  void initialiseInverseSHComputation ( proshade_unsign shBand, double*& sigR, double*& sigI, double*& rcoeffs, double*& icoeffs, double*& weights, double*& workspace,
62  fftw_plan& idctPlan, fftw_plan& ifftPlan );
63 }
64 
65 #endif
ProSHADE_internal_data::ProSHADE_data
This class contains all inputed and derived data for a single structure.
Definition: ProSHADE_data.hpp:49
ProSHADE_internal_overlay::computeBeforeAfterZeroCounts
void computeBeforeAfterZeroCounts(proshade_unsign *addXPre, proshade_unsign *addYPre, proshade_unsign *addZPre, proshade_unsign *addXPost, proshade_unsign *addYPost, proshade_unsign *addZPost, proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim, proshade_unsign xDimIndices, proshade_unsign yDimIndices, proshade_unsign zDimIndices)
This function finds the number of zeroes to be added after and before the structure along each dimens...
Definition: ProSHADE_overlay.cpp:657
ProSHADE_symmetry.hpp
This header file declares all the functions required for symmetry detection and construction.
ProSHADE_internal_overlay::getOptimalRotation
void getOptimalRotation(ProSHADE_settings *settings, ProSHADE_internal_data::ProSHADE_data *staticStructure, ProSHADE_internal_data::ProSHADE_data *movingStructure, proshade_double *eulA, proshade_double *eulB, proshade_double *eulG)
This function finds the optimal rotation between two structures as described by the settings object.
Definition: ProSHADE_overlay.cpp:154
ProSHADE_internal_overlay::freeTranslationFunctionMemory
void freeTranslationFunctionMemory(fftw_complex *&tmpIn1, fftw_complex *&tmpOut1, fftw_complex *&tmpIn2, fftw_complex *&tmpOut2, fftw_complex *&resOut, fftw_plan &forwardFourierObj1, fftw_plan &forwardFourierObj2, fftw_plan &inverseFourierCombo)
This function releases the memory for the Fourier transforms required for translation function comput...
Definition: ProSHADE_overlay.cpp:471
ProSHADE_internal_overlay::getOptimalTranslation
void getOptimalTranslation(ProSHADE_settings *settings, ProSHADE_internal_data::ProSHADE_data *staticStructure, ProSHADE_internal_data::ProSHADE_data *movingStructure, proshade_double *trsX, proshade_double *trsY, proshade_double *trsZ, proshade_double eulA, proshade_double eulB, proshade_double eulG)
This function finds the optimal translation between two structures as described by the settings objec...
Definition: ProSHADE_overlay.cpp:203
ProSHADE_internal_overlay
This namespace contains all the functions required for map overlays.
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:89
ProSHADE_internal_overlay::computeAngularThreshold
void computeAngularThreshold(std::vector< proshade_double > *lonCO, std::vector< proshade_double > *latCO, proshade_unsign angRes)
This function computes the angular thresholds for longitude and lattitude angles.
Definition: ProSHADE_overlay.cpp:1018
ProSHADE_internal_overlay::initialiseInverseSHComputation
void initialiseInverseSHComputation(proshade_unsign shBand, double *&sigR, double *&sigI, double *&rcoeffs, double *&icoeffs, double *&weights, double *&workspace, fftw_plan &idctPlan, fftw_plan &ifftPlan)
This function initialises internal variables for inverse Spherical Harmonics computation.
Definition: ProSHADE_overlay.cpp:902
ProSHADE_internal_overlay::findHighestValueInMap
void findHighestValueInMap(fftw_complex *resIn, proshade_unsign xD, proshade_unsign yD, proshade_unsign zD, proshade_double *trsX, proshade_double *trsY, proshade_double *trsZ, proshade_double *mapPeak)
This function simply finds the highest value in fftw_complex map and returns its position and value.
Definition: ProSHADE_overlay.cpp:558
ProSHADE_internal_overlay::allocateTranslationFunctionMemory
void allocateTranslationFunctionMemory(fftw_complex *&tmpIn1, fftw_complex *&tmpOut1, fftw_complex *&tmpIn2, fftw_complex *&tmpOut2, fftw_complex *&resIn, fftw_complex *&resOut, fftw_plan &forwardFourierObj1, fftw_plan &forwardFourierObj2, fftw_plan &inverseFourierCombo, proshade_unsign xD, proshade_unsign yD, proshade_unsign zD)
This function allocates the memory for the Fourier transforms required for translation function compu...
Definition: ProSHADE_overlay.cpp:432
ProSHADE_internal_overlay::paddMapWithZeroes
void paddMapWithZeroes(proshade_double *oldMap, proshade_double *&newMap, proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim, proshade_unsign xDimIndices, proshade_unsign yDimIndices, proshade_unsign zDimIndices, proshade_unsign addXPre, proshade_unsign addYPre, proshade_unsign addZPre)
This function adds zeroes before and after the central map and copies the central map values into a n...
Definition: ProSHADE_overlay.cpp:686
ProSHADE_internal_overlay::combineFourierForTranslation
void combineFourierForTranslation(fftw_complex *tmpOut1, fftw_complex *tmpOut2, fftw_complex *&resOut, proshade_unsign xD, proshade_unsign yD, proshade_unsign zD)
This function combines Fourier coefficients of two structures in a way, so that inverse Fourier of th...
Definition: ProSHADE_overlay.cpp:497