ProSHADE  0.7.6.6 (JUL 2022)
Protein Shape Detection
ProSHADE_sphericalHarmonics.cpp
Go to the documentation of this file.
1 
22 //==================================================== ProSHADE
24 
39 void ProSHADE_internal_sphericalHarmonics::allocateComputationMemory ( proshade_unsign band, proshade_double*& inputReal, proshade_double*& inputImag, proshade_double*& outputReal, proshade_double*& outputImag, proshade_double*& shWeights, proshade_double*& tableSpaceHelper, fftw_complex*& workspace )
40 {
41  //================================================ Initialise local variables
42  proshade_unsign oneDimmension = 2 * band;
43 
44  //================================================ Allocate Input Memory
45  inputReal = new proshade_double [oneDimmension * oneDimmension];
46  inputImag = new proshade_double [oneDimmension * oneDimmension];
47 
48  //================================================ Allocate Output Memory
49  outputReal = new proshade_double [oneDimmension * oneDimmension];
50  outputImag = new proshade_double [oneDimmension * oneDimmension];
51 
52  //================================================ Allocate Working Memory
53  shWeights = new proshade_double [band * 4];
54  workspace = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * ( 8 * band * band ) + ( 10 * band ) ) );
55 
56  //================================================ Allocate table
57  tableSpaceHelper = new proshade_double [static_cast<proshade_unsign> ( Reduced_Naive_TableSize ( static_cast< int > ( band ), static_cast< int > ( band ) ) +
58  Reduced_SpharmonicTableSize ( static_cast< int > ( band ), static_cast< int > ( band ) ) )];
59 
60  //================================================ Check memory allocation success
61  ProSHADE_internal_misc::checkMemoryAllocation ( inputReal, __FILE__, __LINE__, __func__ );
62  ProSHADE_internal_misc::checkMemoryAllocation ( inputImag, __FILE__, __LINE__, __func__ );
63  ProSHADE_internal_misc::checkMemoryAllocation ( outputReal, __FILE__, __LINE__, __func__ );
64  ProSHADE_internal_misc::checkMemoryAllocation ( outputImag, __FILE__, __LINE__, __func__ );
65  ProSHADE_internal_misc::checkMemoryAllocation ( shWeights, __FILE__, __LINE__, __func__ );
66  ProSHADE_internal_misc::checkMemoryAllocation ( tableSpaceHelper, __FILE__, __LINE__, __func__ );
67  ProSHADE_internal_misc::checkMemoryAllocation ( workspace, __FILE__, __LINE__, __func__ );
68 
69  //================================================ Fill arrays with zeroes
70  for ( size_t iter = 0; iter < static_cast< size_t > ( oneDimmension * oneDimmension ); iter++ )
71  {
72  inputReal[iter] = 0.0;
73  inputImag[iter] = 0.0;
74  outputReal[iter] = 0.0;
75  outputImag[iter] = 0.0;
76  }
77 
78  //================================================ Done
79  return ;
80 
81 }
82 
95 void ProSHADE_internal_sphericalHarmonics::placeWithinWorkspacePointers ( fftw_complex*& workspace, proshade_unsign oDim, proshade_double*& rres, proshade_double*& ires, proshade_double*& fltres, proshade_double*& scratchpad )
96 {
97  //================================================ Place pointers as required by SOFT2.0
98  rres = reinterpret_cast<proshade_double*> ( workspace );
99  ires = rres + ( oDim * oDim );
100  fltres = ires + ( oDim * oDim );
101  scratchpad = fltres + ( oDim / 2 );
102 
103  //================================================ Done
104  return ;
105 
106 }
107 
122 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 )
123 {
124  //================================================ Initialize fft plan along phi angles
125  fftw_iodim dims[1];
126  fftw_iodim howmany_dims[1];
127 
128  int rank = 1;
129  int howmany_rank = 1;
130 
131  dims[0].n = static_cast<int> ( band * 2 );
132  dims[0].is = 1;
133  dims[0].os = static_cast<int> ( band * 2 );
134 
135  howmany_dims[0].n = static_cast<int> ( band * 2 );
136  howmany_dims[0].is = static_cast<int> ( band * 2 );
137  howmany_dims[0].os = 1;
138 
139  //================================================ Plan fft transform
140  fftPlan = fftw_plan_guru_split_dft ( rank,
141  dims,
142  howmany_rank,
143  howmany_dims,
144  inputReal,
145  inputImag,
146  rres,
147  ires,
148  FFTW_ESTIMATE );
149 
150  //================================================ Initialize dct plan for SHT
151  dctPlan = fftw_plan_r2r_1d ( static_cast<int> ( band * 2 ),
152  scratchpad,
153  scratchpad + static_cast<int> ( band * 2 ),
154  FFTW_REDFT10,
155  FFTW_ESTIMATE ) ;
156 
157  //================================================ Done
158  return ;
159 
160 }
161 
162 
179 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 )
180 {
181  //================================================ Release all memory related to SH
182  delete[] inputReal;
183  delete[] inputImag;
184  delete[] outputReal;
185  delete[] outputImag;
186  delete[] tableSpaceHelper;
187  delete[] shWeights;
188  fftw_free ( workspace );
189 
190  //================================================ Set pointers to NULL
191  tableSpaceHelper = nullptr;
192  tableSpace = nullptr;
193  shWeights = nullptr;
194  workspace = nullptr;
195 
196  //================================================ Delete fftw plans
197  fftw_destroy_plan ( dctPlan );
198  fftw_destroy_plan ( fftPlan );
199 
200  //================================================ Done
201  return ;
202 
203 }
204 
226 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 )
227 {
228  //================================================ Initialise local variables
229  proshade_unsign oneDim = band * 2;
230 
231  //================================================ Allocate memory for local pointers
232  allocateComputationMemory ( band, inputReal, inputImag, outputReal, outputImag, shWeights, tableSpaceHelper, workspace );
233 
234  //================================================ Within workspace pointers
235  placeWithinWorkspacePointers ( workspace, oneDim, rres, ires, fltres, scratchpad );
236 
237  //================================================ Generate Seminaive and naive tables for Legendre Polynomials
238  tableSpace = SemiNaive_Naive_Pml_Table ( static_cast< int > ( band ), static_cast< int > ( band ), tableSpaceHelper, reinterpret_cast<double*> ( workspace ) );
239 
240  //================================================ Make weights for spherical transform
241  makeweights ( static_cast< int > ( band ), shWeights );
242 
243  //================================================ Initialize FFTW Plans
244  initialiseFFTWPlans ( band, fftPlan, dctPlan, inputReal, inputImag, rres, ires, scratchpad );
245 
246  //================================================ Done
247  return ;
248 
249 }
250 
266 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 )
267 {
268  //================================================ Load mapped data to decomposition array
269  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( oneDim * oneDim ); iter++ )
270  {
271  inputReal[iter] = mappedData[iter];
272  inputImag[iter] = 0.0;
273  }
274 
275  //================================================ Execute fft plan along phi
276  fftw_execute_split_dft ( fftPlan, inputReal, inputImag, rres, ires ) ;
277 
278  //================================================ Normalize
279  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( oneDim * oneDim ); iter++ )
280  {
281  rres[iter] *= normCoeff;
282  ires[iter] *= normCoeff;
283  }
284 
285  //================================================ Done
286  return ;
287 
288 }
289 
309 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 )
310 {
311  //================================================ Calculate the coefficients for each band
312  rdataptr = outputReal;
313  idataptr = outputImag;
314  for ( proshade_unsign bandIter = 0; bandIter < band; bandIter++ )
315  {
316  //============================================ Real part calculation
317  SemiNaiveReduced ( rres + ( bandIter * ( band * 2 ) ),
318  static_cast< int > ( band ),
319  static_cast< int > ( bandIter ),
320  fltres,
321  scratchpad,
322  tablePml[bandIter],
323  shWeights,
324  &dctPlan);
325 
326  //============================================ Save the real results to temporary holder
327  memcpy ( rdataptr, fltres, sizeof(proshade_double) * ( band - bandIter ) );
328  rdataptr += band - bandIter;
329 
330  //============================================ Imaginary part calculation
331  SemiNaiveReduced ( ires + ( bandIter * ( band * 2 ) ),
332  static_cast< int > ( band ),
333  static_cast< int > ( bandIter ),
334  fltres,
335  scratchpad,
336  tablePml[bandIter],
337  shWeights,
338  &dctPlan );
339 
340  //============================================ Save the imaginary results
341  memcpy ( idataptr, fltres, sizeof(proshade_double) * ( band - bandIter ) );
342  idataptr += band - bandIter;
343  }
344 
345  //================================================ DONE
346  return ;
347 
348 }
349 
361 void ProSHADE_internal_sphericalHarmonics::applyCondonShortleyPhase ( proshade_unsign band, proshade_double* outputReal, proshade_double* outputImag, proshade_complex*& shArray )
362 {
363  //================================================ Copy the results into the final holder
364  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( (band * 2) * (band * 2) ); iter++ )
365  {
366  shArray[iter][0] = outputReal[iter];
367  shArray[iter][1] = outputImag[iter];
368  }
369 
370  //================================================ Apply the Condon-Shortley phase sign
371  proshade_double powerOne = 1.0;
372  proshade_unsign hlp1 = 0;
373  proshade_unsign hlp2 = 0;
374  for ( proshade_signed order = 1; order < static_cast<proshade_signed> ( band ); order++)
375  {
376  powerOne *= -1.0;
377  for ( proshade_signed bandIter = order; bandIter < static_cast<proshade_signed> ( band ); bandIter++)
378  {
379  hlp1 = static_cast< proshade_unsign > ( seanindex ( static_cast< int > ( order ), static_cast< int > ( bandIter ), static_cast< int > ( band ) ) );
380  hlp2 = static_cast< proshade_unsign > ( seanindex ( static_cast< int > ( -order ), static_cast< int > ( bandIter ), static_cast< int > ( band ) ) );
381 
382  shArray[hlp2][0] = powerOne * static_cast<proshade_double> ( outputReal[hlp1] );
383  shArray[hlp2][1] = -powerOne * static_cast<proshade_double> ( outputImag[hlp1] );
384  }
385  }
386 
387  //================================================ DONE
388  return ;
389 
390 }
391 
403 void ProSHADE_internal_sphericalHarmonics::computeSphericalHarmonics ( proshade_unsign band, proshade_double* sphereMappedData, proshade_complex*& shArray )
404 {
405  //================================================ Initialise local variables
406  proshade_double *inputReal = nullptr, *inputImag = nullptr, *outputReal = nullptr, *outputImag = nullptr;
407  double *shWeights = nullptr, *tableSpaceHelper = nullptr;
408  double** tablePml = nullptr;
409  fftw_complex* workspace = nullptr;
410  proshade_unsign oneDim = static_cast<proshade_unsign> ( band * 2 );
411  proshade_double normCoeff = ( 1.0 / ( static_cast<proshade_double> ( band * 2 ) ) ) * sqrt( 2.0 * M_PI );
412 
413  //================================================ Set output to zeroes (so that all unfilled data are not random)
414  for ( proshade_unsign i = 0; i < ( 2 * band * 2 * band); i++ )
415  {
416  shArray[i][0] = 0.0;
417  shArray[i][1] = 0.0;
418  }
419 
420  //================================================ Within workspace pointers
421  proshade_double *rres = nullptr, *ires = nullptr, *fltres = nullptr, *scratchpad = nullptr, *rdataptr = nullptr, *idataptr = nullptr;
422 
423  //================================================ FFTW Plans
424  fftw_plan fftPlan = nullptr;
425  fftw_plan dctPlan = nullptr;
426 
427  //================================================ Initialise all memory
428  initialiseAllMemory ( band, inputReal, inputImag, outputReal, outputImag, shWeights, tablePml, tableSpaceHelper, workspace,
429  rres, ires, fltres, scratchpad, fftPlan, dctPlan );
430 
431  //================================================ Do the initial discrete split transform
432  initialSplitDiscreteTransform ( oneDim, inputReal, inputImag, rres, ires, sphereMappedData, fftPlan, normCoeff );
433 
434  //================================================ Complete the spherical harmonics transform
435  computeSphericalTransformCoeffs ( band, rdataptr, idataptr, outputReal, outputImag, rres, ires, fltres, scratchpad, tablePml, shWeights, dctPlan );
436 
437  //================================================ Apply the Condon-Shortley phase and save result to the final array
438  applyCondonShortleyPhase ( band, outputReal, outputImag, shArray );
439 
440  //================================================ Free memory
441  releaseSphericalMemory ( inputReal, inputImag, outputReal, outputImag, tableSpaceHelper, tablePml, shWeights, workspace, fftPlan, dctPlan );
442 
443  //================================================ Done
444  return ;
445 
446 }
ProSHADE_internal_sphericalHarmonics::initialiseAllMemory
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.
Definition: ProSHADE_sphericalHarmonics.cpp:226
ProSHADE_internal_sphericalHarmonics::allocateComputationMemory
void allocateComputationMemory(proshade_unsign band, proshade_double *&inputReal, proshade_double *&inputImag, proshade_double *&outputReal, proshade_double *&outputImag, double *&shWeights, double *&tableSpaceHelper, fftw_complex *&workspace)
This function determines the integration order for the between spheres integration.
Definition: ProSHADE_sphericalHarmonics.cpp:39
ProSHADE_internal_sphericalHarmonics::releaseSphericalMemory
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)
This function computes the spherical harmonics of a aingle shell, saving them in supplied pointer.
Definition: ProSHADE_sphericalHarmonics.cpp:179
ProSHADE_internal_sphericalHarmonics::initialSplitDiscreteTransform
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.
Definition: ProSHADE_sphericalHarmonics.cpp:266
ProSHADE_internal_sphericalHarmonics::computeSphericalTransformCoeffs
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 dec...
Definition: ProSHADE_sphericalHarmonics.cpp:309
ProSHADE_internal_misc::checkMemoryAllocation
void checkMemoryAllocation(chVar checkVar, std::string fileP, unsigned int lineP, std::string funcP, std::string infoP="This error may occurs when ProSHADE requests memory to be\n : allocated to it and this operation fails. This could\n : happen when not enough memory is available, either due to\n : other processes using a lot of memory, or when the machine\n : does not have sufficient memory available. Re-run to see\n : if this problem persists.")
Checks if memory was allocated properly.
Definition: ProSHADE_misc.hpp:73
ProSHADE_internal_sphericalHarmonics::placeWithinWorkspacePointers
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.
Definition: ProSHADE_sphericalHarmonics.cpp:95
ProSHADE_internal_sphericalHarmonics::applyCondonShortleyPhase
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.
Definition: ProSHADE_sphericalHarmonics.cpp:361
ProSHADE_internal_sphericalHarmonics::initialiseFFTWPlans
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.
Definition: ProSHADE_sphericalHarmonics.cpp:122
ProSHADE_internal_sphericalHarmonics::computeSphericalHarmonics
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.
Definition: ProSHADE_sphericalHarmonics.cpp:403
ProSHADE_sphericalHarmonics.hpp
This header file declares the functions required to compute the spherical harmonics decomposition for...