![]() |
ProSHADE
0.7.5.3 (FEB 2021)
Protein Shape Detection
|
This namespace contains the internal functions for computing spherical harmonics and their related computations. More...
Functions | |
void | allocateComputationMemory (proshade_unsign band, proshade_double *&inputReal, proshade_double *&inputImag, proshade_double *&outputReal, proshade_double *&outputImag, double *&shWeights, double **&tableSpace, double *&tableSpaceHelper, fftw_complex *&workspace) |
This function determines the integration order for the between spheres integration. More... | |
void | placeWithinWorkspacePointers (fftw_complex *&workspace, proshade_unsign oDim, proshade_double *&rres, proshade_double *&ires, proshade_double *&fltres, proshade_double *&scratchpad) |
This function takes the workspace pointer and correctly places the other internal pointers. More... | |
void | initialiseFFTWPlans (proshade_unsign band, fftw_plan &fftPlan, fftw_plan &dctPlan, proshade_double *&inputReal, proshade_double *&inputImag, proshade_double *&rres, proshade_double *&ires, proshade_double *&scratchpad) |
This function initialises the FFTW plans. More... | |
void | releaseSphericalMemory (proshade_double *&inputReal, proshade_double *&inputImag, proshade_double *&outputReal, proshade_double *&outputImag, double *&tableSpaceHelper, double **&tableSpace, double *&shWeights, fftw_complex *&workspace, fftw_plan &fftPlan, fftw_plan &dctPlan, proshade_unsign band) |
This function computes the spherical harmonics of a aingle shell, saving them in supplied pointer. More... | |
void | initialiseAllMemory (proshade_unsign band, proshade_double *&inputReal, proshade_double *&inputImag, proshade_double *&outputReal, proshade_double *&outputImag, double *&shWeights, double **&tableSpace, double *&tableSpaceHelper, fftw_complex *&workspace, proshade_double *&rres, proshade_double *&ires, proshade_double *&fltres, proshade_double *&scratchpad, fftw_plan &fftPlan, fftw_plan &dctPlan) |
This function initialises all the memory required for spherical harmonics computation. More... | |
void | initialSplitDiscreteTransform (proshade_unsign oneDim, proshade_double *&inputReal, proshade_double *&inputImag, proshade_double *&rres, proshade_double *&ires, proshade_double *mappedData, fftw_plan &fftPlan, proshade_double normCoeff) |
This function computes the spherical harmonics of a aingle shell, saving them in supplied pointer. More... | |
void | computeSphericalTransformCoeffs (proshade_unsign band, proshade_double *&rdataptr, proshade_double *&idataptr, proshade_double *&outputReal, proshade_double *&outputImag, proshade_double *&rres, proshade_double *&ires, proshade_double *&fltres, proshade_double *&scratchpad, double **&tablePml, double *&shWeights, fftw_plan &dctPlan) |
This function takes the split discrete transform and proceeds to complete the spherical harmonics decomposition. More... | |
void | applyCondonShortleyPhase (proshade_unsign band, proshade_double *outputReal, proshade_double *outputImag, proshade_complex *&shArray) |
This is the final step in computing the full spherical harmonics decomposition of the input data. More... | |
void | computeSphericalHarmonics (proshade_unsign band, proshade_double *sphereMappedData, proshade_complex *&shArray) |
This function computes the spherical harmonics of a aingle shell, saving them in supplied pointer. More... | |
This namespace contains the internal functions for computing spherical harmonics and their related computations.
The ProSHADE_internal_sphericalHarmonics namespace contains helper functions for spherical harmonics computation for each sphere as created in the sphere object as well as the further processing computations. None of these functions should be used directly be the user.
void ProSHADE_internal_sphericalHarmonics::allocateComputationMemory | ( | proshade_unsign | band, |
proshade_double *& | inputReal, | ||
proshade_double *& | inputImag, | ||
proshade_double *& | outputReal, | ||
proshade_double *& | outputImag, | ||
double *& | shWeights, | ||
double **& | tableSpace, | ||
double *& | tableSpaceHelper, | ||
fftw_complex *& | workspace | ||
) |
This function determines the integration order for the between spheres integration.
This function simply takes all pointer variables required for the spherical harmonics computation and allocates the required amount of memory for them. It also does the memory checks in case memory allocation fails.
[in] | band | The bandwidth to which the computation will be done. |
[in] | inputReal | The real part of the input will be copied here. |
[in] | inputImag | The immaginary part of the input will be copied here. |
[in] | outputReal | The real part of the output will be saved here. |
[in] | outputImag | The immaginary part of the output will be saved here. |
[in] | shWeights | The weights for spherical harmonics computation will be stored here. |
[in] | tableSpace | Space where the legendre polynomials table will be stored. |
[in] | tableSpaceHelper | This space is required by SOFT for pre-computing values into this table. |
[in] | workspace | The space where multiple minor results are saved by SOFT. |
Definition at line 40 of file ProSHADE_sphericalHarmonics.cpp.
void ProSHADE_internal_sphericalHarmonics::applyCondonShortleyPhase | ( | proshade_unsign | band, |
proshade_double * | outputReal, | ||
proshade_double * | outputImag, | ||
proshade_complex *& | shArray | ||
) |
This is the final step in computing the full spherical harmonics decomposition of the input data.
This is the final function that is needed for the complete spherical harmonics decomposition. It applied the Condon-Shortley phase as well as computing the negative orders, then saving the output into the final results array for further processing.
[in] | band | The bandwidth to which the computation will be done. |
[in] | outputReal | The real results of the complete transform as done by the initialSplitDiscreteTransform() and computeSphericalTransformCoeffs() functions. |
[in] | outputImag | The imaginary results of the complete transform as done by the initialSplitDiscreteTransform() and computeSphericalTransformCoeffs() functions. |
[in] | shArray | An array of complex numbers to which the results of the spherical harmonics decomposition are to be saved. |
Definition at line 353 of file ProSHADE_sphericalHarmonics.cpp.
void ProSHADE_internal_sphericalHarmonics::computeSphericalHarmonics | ( | proshade_unsign | band, |
proshade_double * | sphereMappedData, | ||
proshade_complex *& | shArray | ||
) |
This function computes the spherical harmonics of a aingle shell, saving them in supplied pointer.
This function does all the spherical harmonics computations for a single shell, including the memory allocation and releasing and the FFTW transforms. Because the shells can have different resolutions, the memory management is left until here, but possible speed-up could be gained from having the same resolution on all shells as in the older ProSHADE versions.
[in] | band | The bandwidth to which the computation will be done. |
[in] | sphereMappedData | An array of doubles containing the mapped data onto a sphere for the sphere to be decomposed. |
[in] | shArray | An array of complex numbers to which the results of the spherical harmonics decomposition are to be saved. |
Definition at line 395 of file ProSHADE_sphericalHarmonics.cpp.
void ProSHADE_internal_sphericalHarmonics::computeSphericalTransformCoeffs | ( | proshade_unsign | band, |
proshade_double *& | rdataptr, | ||
proshade_double *& | idataptr, | ||
proshade_double *& | outputReal, | ||
proshade_double *& | outputImag, | ||
proshade_double *& | rres, | ||
proshade_double *& | ires, | ||
proshade_double *& | fltres, | ||
proshade_double *& | scratchpad, | ||
double **& | tablePml, | ||
double *& | shWeights, | ||
fftw_plan & | dctPlan | ||
) |
This function takes the split discrete transform and proceeds to complete the spherical harmonics decomposition.
This function takes the results of the initial split discrete transform and proceeds to compute the spherical harmonics coefficients for all applicable bands using all the pre-computed valus (i.e. the Legendre polynomials table and the weights). To do this, the SOFT2.0 function is called and for more details, see SOFT2.0.
[in] | band | The bandwidth to which the computation will be done. |
[in] | rdataptr | Pointer to be used as iterative locator in the large results array for real results. |
[in] | idataptr | Pointer to be used as iterative locator in the large results array for imaginary results. |
[in] | outputReal | An array for storing the real results of the completed transform. |
[in] | outputImag | An array for storing the imaginary results of the completed transform. |
[in] | rres | Array containing the real results of the initial split discrete transform. |
[in] | ires | Array containing the imaginary results of the initial split discrete transform. |
[in] | fltres | Helper array for transform computations. |
[in] | scratchpad | Array for keeping temporary results of the transform computations. |
[in] | tablePml | Pre-computed array of the Legendre polynomials as done by SOFT2.0. |
[in] | shWeights | The weights for the spherical harmonics. |
[in] | dctPlan | The FFTW plan for the final spherical harmonics transform. |
Definition at line 301 of file ProSHADE_sphericalHarmonics.cpp.
void ProSHADE_internal_sphericalHarmonics::initialiseAllMemory | ( | proshade_unsign | band, |
proshade_double *& | inputReal, | ||
proshade_double *& | inputImag, | ||
proshade_double *& | outputReal, | ||
proshade_double *& | outputImag, | ||
double *& | shWeights, | ||
double **& | tableSpace, | ||
double *& | tableSpaceHelper, | ||
fftw_complex *& | workspace, | ||
proshade_double *& | rres, | ||
proshade_double *& | ires, | ||
proshade_double *& | fltres, | ||
proshade_double *& | scratchpad, | ||
fftw_plan & | fftPlan, | ||
fftw_plan & | dctPlan | ||
) |
This function initialises all the memory required for spherical harmonics computation.
This function takes on all the memory allocation and filling in all data required later for the spherical harmonics computation by the SOFT2.0 library.
[in] | band | The bandwidth to which the computation will be done. |
[in] | inputReal | The real part of the input will be copied here. |
[in] | inputImag | The immaginary part of the input will be copied here. |
[in] | outputReal | The real part of the output will be saved here. |
[in] | outputImag | The immaginary part of the output will be saved here. |
[in] | shWeights | The weights for spherical harmonics computation will be stored here. |
[in] | tableSpace | This space is required by SOFT for pre-computing values into this table. |
[in] | tableSpaceHelper | This space is required by SOFT for proper computation of the table space. |
[in] | workspace | The space where multiple minor results are saved by SOFT2.0. |
[in] | rres | Pointer to where the real part of the results will be temporarily saved. |
[in] | ires | Pointer to where the imaginary part of the results will be temporarily saved. |
[in] | fltres | Pointer to where the temporary bandwise results will be saved. |
[in] | scratchpad | Pointer to extra space which is internally used (but not allocated) by SOFT2.0. |
[in] | fftPlan | pointer to the variable where the Fourier transform should be set. |
[in] | dctPlan | pointer to the variable where the 1D r2r Fourier transform should be set. |
Definition at line 218 of file ProSHADE_sphericalHarmonics.cpp.
void ProSHADE_internal_sphericalHarmonics::initialiseFFTWPlans | ( | proshade_unsign | band, |
fftw_plan & | fftPlan, | ||
fftw_plan & | dctPlan, | ||
proshade_double *& | inputReal, | ||
proshade_double *& | inputImag, | ||
proshade_double *& | rres, | ||
proshade_double *& | ires, | ||
proshade_double *& | scratchpad | ||
) |
This function initialises the FFTW plans.
This function initialises the FFTW plans for the spherical harmonics computations as required by the SOFT2.0 library.
[in] | band | The bandwidth to which the computation will be done. |
[in] | fftPlan | pointer to the variable where the Fourier transform should be set. |
[in] | dctPlan | pointer to the variable where the 1D r2r Fourier transform should be set. |
[in] | inputReal | pointer to the array containing (or which will contain) the input real values. |
[in] | inputImag | pointer to the array containing (or which will contain) the input imaginary values. |
[in] | rres | pointer to the array where the real values of result should be saved. |
[in] | ires | pointer to the array where the imaginary values of result should be saved. |
[in] | scratchpad | pointer to the array where temporary results will be saved. |
Definition at line 114 of file ProSHADE_sphericalHarmonics.cpp.
void ProSHADE_internal_sphericalHarmonics::initialSplitDiscreteTransform | ( | proshade_unsign | oneDim, |
proshade_double *& | inputReal, | ||
proshade_double *& | inputImag, | ||
proshade_double *& | rres, | ||
proshade_double *& | ires, | ||
proshade_double * | mappedData, | ||
fftw_plan & | fftPlan, | ||
proshade_double | normCoeff | ||
) |
This function computes the spherical harmonics of a aingle shell, saving them in supplied pointer.
This function takes the already initialised and prepared pointers and values and proceeds to load the data into the proper places and compute the split discrete Fourier transform, thus preparing for the spherical transform to be done.
[in] | oneDim | This is the size of any dimension of the transform (2 * bandwidth). |
[in] | inputReal | Pointer to array which should be subjected to the transform (real part). |
[in] | inputReal | Pointer to array which should be subjected to the transform (imaginary part). |
[in] | rres | Pointer to array where the transform results will be saved (real part). |
[in] | ires | Pointer to array where the transform results will be saved (imaginary part). |
[in] | mappedData | Pointer to the data which should be decomposed. |
[in] | fftPlan | The prepared plan which states how the transform will be done as set by the initialiseFFTWPlans() function. |
[in] | normCoeff | The transform normalisation factor. |
Definition at line 258 of file ProSHADE_sphericalHarmonics.cpp.
void ProSHADE_internal_sphericalHarmonics::placeWithinWorkspacePointers | ( | fftw_complex *& | workspace, |
proshade_unsign | oDim, | ||
proshade_double *& | rres, | ||
proshade_double *& | ires, | ||
proshade_double *& | fltres, | ||
proshade_double *& | scratchpad | ||
) |
This function takes the workspace pointer and correctly places the other internal pointers.
This is a simple helper function, which places the internal workspace pointers to the correct addresses as required by the SOFT2.0 dependency.
[in] | workspace | Pointer to the allocated workspace, within which the other pointers will be placed. |
[in] | oDim | The size of the single transform dimension (twice the transform bandwidth). |
[in] | rres | Pointer to where the real part of the results will be temporarily saved. |
[in] | ires | Pointer to where the imaginary part of the results will be temporarily saved. |
[in] | fltres | Pointer to where the temporary bandwise results will be saved. |
[in] | scratchpad | Pointer to extra space which is internally used (but not allocated) by SOFT2.0. |
Definition at line 87 of file ProSHADE_sphericalHarmonics.cpp.
void ProSHADE_internal_sphericalHarmonics::releaseSphericalMemory | ( | proshade_double *& | inputReal, |
proshade_double *& | inputImag, | ||
proshade_double *& | outputReal, | ||
proshade_double *& | outputImag, | ||
double *& | tableSpaceHelper, | ||
double **& | tableSpace, | ||
double *& | shWeights, | ||
fftw_complex *& | workspace, | ||
fftw_plan & | fftPlan, | ||
fftw_plan & | dctPlan, | ||
proshade_unsign | band | ||
) |
This function computes the spherical harmonics of a aingle shell, saving them in supplied pointer.
This function takes all the pointers required for the spherical harmonics computation and releases all the memory.
[in] | inputReal | pointer to the array that contained the input real values to be freed. |
[in] | inputImag | pointer to the array that contained the input imaginary values to be freed. |
[in] | outputReal | pointer to the array that contained the output real values to be freed. |
[in] | outputImag | pointer to the array that contained the output imaginary values to be freed. |
[in] | tableSpaceHelper | pointer to the helper array for Legendre polynomials values to be freed. |
[in] | tableSpace | pointer to the array of Legendre polynomials values to be freed. |
[in] | shWeights | pointer to the array spherical harmonics weighhts to be freed. |
[in] | workspace | pointer to the array for miscellaneous temporary results to be freed. |
[in] | fftPlan | pointer to the variable where the Fourier transform was done to be freed. |
[in] | dctPlan | pointer to the variable where the 1D r2r Fourier transform was done to be freed. |
Definition at line 171 of file ProSHADE_sphericalHarmonics.cpp.