ProSHADE  0.6.6 (DEC 2018)
Protein Shape Descriptors and Symmetry Detection
ProSHADE_spharmonickit.cpp
Go to the documentation of this file.
1 
20 //============================================ FFTW3 + SOFT
21 #ifdef __cplusplus
22 extern "C" {
23 #endif
24 #include <fftw3.h>
25 #include <wrap_fftw.h>
26 #include <makeweights.h>
27 #include <s2_primitive.h>
28 #include <s2_cospmls.h>
29 #include <s2_legendreTransforms.h>
30 #include <s2_semi_fly.h>
31 #include <rotate_so3_utils.h>
32 #include <utils_so3.h>
33 #include <soft_fftw.h>
34 #include <rotate_so3_fftw.h>
35 #ifdef __cplusplus
36 }
37 #endif
38 
39 //============================================ RVAPI
40 #include <rvapi_interface.h>
41 
42 //============================================ ProSHADE
43 #include "ProSHADE.h"
44 #include "ProSHADE_internal.h"
45 
61  double theta,
62  double phi,
63  double shellSz,
64  unsigned int manualShells,
65  bool keepInMemory,
66  bool rotDefaults )
67 {
68  //======================================== Sanity check
69  if ( !this->_phaseRemoved )
70  {
71  std::cerr << "!!! ProSHADE ERROR !!! Error in file " << this->_inputFileName << " !!! Attempted to map data onto shell before computing the phaseless data. Call the removePhaseFromMap function BEFORE the mapPhaselessToSphere function." << std::endl;
72 
73  if ( settings->htmlReport )
74  {
75  std::stringstream hlpSS;
76  hlpSS << "<font color=\"red\">" << "Attempted to map data onto shell before computing the phaseless data. This looks like an internal bug, please report this case." << "</font>";
77  rvapi_set_text ( hlpSS.str().c_str(),
78  "ProgressSection",
79  settings->htmlReportLineProgress,
80  1,
81  1,
82  1 );
83  settings->htmlReportLineProgress += 1;
84  rvapi_flush ( );
85  }
86 
87  exit ( -1 );
88  }
89 
90  //======================================== Initialise internal values
91  this->_thetaAngle = theta;
92  this->_phiAngle = phi;
93  this->_xSamplingRate = this->_xRange / static_cast<double> ( this->_maxMapU );
94  this->_ySamplingRate = this->_yRange / static_cast<double> ( this->_maxMapV );
95  this->_zSamplingRate = this->_zRange / static_cast<double> ( this->_maxMapW );
96 
97  //======================================== Initialise local variables
98  double x, y, z;
99  int xBottom, yBottom, zBottom;
100  int xTop, yTop, zTop;
101 
102  std::array<double,4> c000, c001, c010, c011, c100, c101, c110, c111;
103  std::array<double,4> c00, c01, c10, c11;
104  std::array<double,4> c0, c1;
105  double xd, yd, zd;
106 
107  //======================================== How many shells can the map fill?
108  this->_shellSpacing = shellSz;
109 
110  unsigned int shellNoHlp = 0;
111  double maxDist = std::max ( sqrt( pow( ( (((this->_maxMapU+1)/2)-1) * this->_xSamplingRate ), 2.0 ) +
112  pow ( ( (((this->_maxMapV+1)/2)-1) * this->_ySamplingRate ), 2.0 ) ),
113  std::max ( sqrt( pow( ( (((this->_maxMapU+1)/2)-1) * this->_xSamplingRate ), 2.0 ) +
114  pow ( ( (((this->_maxMapW+1)/2)-1) * this->_zSamplingRate ), 2.0 ) ),
115  sqrt( pow( ( (((this->_maxMapW+1)/2)-1) * this->_zSamplingRate ), 2.0 ) +
116  pow ( ( (((this->_maxMapV+1)/2)-1) * this->_ySamplingRate ), 2.0 ) ) ) );
117  for ( unsigned int sh = 0; sh < static_cast<unsigned int> ( this->_shellPlacement.size() ); sh++ )
118  {
119  if ( maxDist > this->_shellPlacement.at(sh) ) { shellNoHlp = sh; }
120  else { break; }
121  }
122  shellNoHlp -= 1;
123 
124  if ( manualShells == 0 ) { this->_noShellsWithData = shellNoHlp; }
125  else { this->_noShellsWithData = manualShells; }
126 
127  if ( this->_noShellsWithData < 2 )
128  {
129  std::cerr << "!!! ProSHADE ERROR !!! Did not find enough shells to which the data could be mapped. Either the structure is rather small and the default values do not work - use the \'sphDistance\' option to adjust, or this is a bug. Terminating..." << std::endl;
130 
131  if ( settings->htmlReport )
132  {
133  std::stringstream hlpSS;
134  hlpSS << "<font color=\"red\">" << "Did not find enough shells to which the data could be mapped. Either the structure is rather small and the default values do not work - use the \'sphDistance\' option to adjust." << "</font>";
135  rvapi_set_text ( hlpSS.str().c_str(),
136  "ProgressSection",
137  settings->htmlReportLineProgress,
138  1,
139  1,
140  1 );
141  settings->htmlReportLineProgress += 1;
142  rvapi_flush ( );
143  }
144 
145  exit ( -1 );
146  }
147 
148  //======================================== Initialise internal values ( and allocate memory )
149  this->_shellMappedData = new double*[this->_noShellsWithData];
150  for ( unsigned int sh = 0; sh < this->_noShellsWithData; sh++ )
151  {
152  this->_shellMappedData[sh] = new double[static_cast<unsigned int> ( this->_thetaAngle * this->_phiAngle )];
153  }
154 
155  //======================================== Find pixelisation cutOffs
156  std::vector<double> lonCO ( this->_thetaAngle + 1 );
157  for ( unsigned int iter = 0; iter <= this->_thetaAngle; iter++ ) { lonCO.at(iter) = static_cast<double> ( iter ) * ( ( static_cast<double> ( M_PI ) * 2.0 ) / static_cast<double> ( this->_thetaAngle ) ) - ( static_cast<double> ( M_PI ) ); }
158  std::vector<double> latCO ( this->_phiAngle + 1 );
159  for ( unsigned int iter = 0; iter <= this->_phiAngle; iter++ ) { latCO.at(iter) = ( static_cast<double> ( iter ) * ( static_cast<double> ( M_PI ) / static_cast<double> ( this->_phiAngle ) ) - ( static_cast<double> ( M_PI ) / 2.0 ) ); }
160 
161  //======================================== Tri-linear interpolation
162  unsigned int posIter;
163 
164  for ( unsigned int sh = 0; sh < this->_noShellsWithData; sh++ )
165  {
166  //==================================== Find pixel centres
167  for ( unsigned int thIt = 0; thIt < this->_thetaAngle; thIt++ )
168  {
169  for ( unsigned int phIt = 0; phIt < this->_phiAngle; phIt++ )
170  {
171  //============================ Get grid point x, y and z
172  x = this->_shellPlacement.at(sh) * cos( ( static_cast<double> ( lonCO.at(thIt) ) + static_cast<double> ( lonCO.at(thIt+1) ) ) / 2.0 ) * cos( ( static_cast<double> ( latCO.at(phIt) ) + static_cast<double> ( latCO.at (phIt+1) ) ) / 2.0 );
173  y = this->_shellPlacement.at(sh) * sin( ( static_cast<double> ( lonCO.at(thIt) ) + static_cast<double> ( lonCO.at(thIt+1) ) ) / 2.0 ) * cos( ( static_cast<double> ( latCO.at(phIt) ) + static_cast<double> ( latCO.at (phIt+1) ) ) / 2.0 );
174  z = this->_shellPlacement.at(sh) * sin( ( static_cast<double> ( latCO.at(phIt) ) + static_cast<double> ( latCO.at(phIt+1) ) ) / 2.0 );
175 
176  //============================ Find 8 closest point around the grid point
177  xBottom = floor ( (x / this->_xSamplingRate) ) + ((this->_maxMapU+1)/2);
178  yBottom = floor ( (y / this->_ySamplingRate) ) + ((this->_maxMapV+1)/2);
179  zBottom = floor ( (z / this->_zSamplingRate) ) + ((this->_maxMapW+1)/2);
180 
181  xTop = xBottom + 1;
182  yTop = yBottom + 1;
183  zTop = zBottom + 1;
184 
185  posIter = zBottom + (this->_maxMapW + 1) * ( yBottom + (this->_maxMapV + 1) * xBottom );
186  if ( posIter > ( (this->_maxMapU+1) * (this->_maxMapV+1) * (this->_maxMapW+1) ) ) { this->_shellMappedData[sh][static_cast<unsigned int> ( phIt * this->_thetaAngle + thIt)] = 0.0; continue; }
187  c000[0] = (this->_densityMapCorCoords[posIter])[0] * this->_xSamplingRate;
188  c000[1] = (this->_densityMapCorCoords[posIter])[1] * this->_ySamplingRate;
189  c000[2] = (this->_densityMapCorCoords[posIter])[2] * this->_zSamplingRate;
190  c000[3] = this->_densityMapCor[posIter];
191 
192  posIter = zTop + (this->_maxMapW + 1) * ( yBottom + (this->_maxMapV + 1) * xBottom );
193  if ( posIter > ( (this->_maxMapU+1) * (this->_maxMapV+1) * (this->_maxMapW+1) ) ) { this->_shellMappedData[sh][static_cast<unsigned int> ( phIt * this->_thetaAngle + thIt)] = 0.0; continue; }
194  c001[0] = (this->_densityMapCorCoords[posIter])[0] * this->_xSamplingRate;
195  c001[1] = (this->_densityMapCorCoords[posIter])[1] * this->_ySamplingRate;
196  c001[2] = (this->_densityMapCorCoords[posIter])[2] * this->_zSamplingRate;
197  c001[3] = this->_densityMapCor[posIter];
198 
199  posIter = zBottom + (this->_maxMapW + 1) * ( yTop + (this->_maxMapV + 1) * xBottom );
200  if ( posIter > ( (this->_maxMapU+1) * (this->_maxMapV+1) * (this->_maxMapW+1) ) ) { this->_shellMappedData[sh][static_cast<unsigned int> ( phIt * this->_thetaAngle + thIt)] = 0.0; continue; }
201  c010[0] = (this->_densityMapCorCoords[posIter])[0] * this->_xSamplingRate;
202  c010[1] = (this->_densityMapCorCoords[posIter])[1] * this->_ySamplingRate;
203  c010[2] = (this->_densityMapCorCoords[posIter])[2] * this->_zSamplingRate;
204  c010[3] = this->_densityMapCor[posIter];
205 
206  posIter = zTop + (this->_maxMapW + 1) * ( yTop + (this->_maxMapV + 1) * xBottom );
207  if ( posIter > ( (this->_maxMapU+1) * (this->_maxMapV+1) * (this->_maxMapW+1) ) ) { this->_shellMappedData[sh][static_cast<unsigned int> ( phIt * this->_thetaAngle + thIt)] = 0.0; continue; }
208  c011[0] = (this->_densityMapCorCoords[posIter])[0] * this->_xSamplingRate;
209  c011[1] = (this->_densityMapCorCoords[posIter])[1] * this->_ySamplingRate;
210  c011[2] = (this->_densityMapCorCoords[posIter])[2] * this->_zSamplingRate;
211  c011[3] = this->_densityMapCor[posIter];
212 
213  posIter = zBottom + (this->_maxMapW + 1) * ( yBottom + (this->_maxMapV + 1) * xTop );
214  if ( posIter > ( (this->_maxMapU+1) * (this->_maxMapV+1) * (this->_maxMapW+1) ) ) { this->_shellMappedData[sh][static_cast<unsigned int> ( phIt * this->_thetaAngle + thIt)] = 0.0; continue; }
215  c100[0] = (this->_densityMapCorCoords[posIter])[0] * this->_xSamplingRate;
216  c100[1] = (this->_densityMapCorCoords[posIter])[1] * this->_ySamplingRate;
217  c100[2] = (this->_densityMapCorCoords[posIter])[2] * this->_zSamplingRate;
218  c100[3] = this->_densityMapCor[posIter];
219 
220  posIter = zTop + (this->_maxMapW + 1) * ( yBottom + (this->_maxMapV + 1) * xTop );
221  if ( posIter > ( (this->_maxMapU+1) * (this->_maxMapV+1) * (this->_maxMapW+1) ) ) { this->_shellMappedData[sh][static_cast<unsigned int> ( phIt * this->_thetaAngle + thIt)] = 0.0; continue; }
222  c101[0] = (this->_densityMapCorCoords[posIter])[0] * this->_xSamplingRate;
223  c101[1] = (this->_densityMapCorCoords[posIter])[1] * this->_ySamplingRate;
224  c101[2] = (this->_densityMapCorCoords[posIter])[2] * this->_zSamplingRate;
225  c101[3] = this->_densityMapCor[posIter];
226 
227  posIter = zBottom + (this->_maxMapW + 1) * ( yTop + (this->_maxMapV + 1) * xTop );
228  if ( posIter > ( (this->_maxMapU+1) * (this->_maxMapV+1) * (this->_maxMapW+1) ) ) { this->_shellMappedData[sh][static_cast<unsigned int> ( phIt * this->_thetaAngle + thIt)] = 0.0; continue; }
229  c110[0] = (this->_densityMapCorCoords[posIter])[0] * this->_xSamplingRate;
230  c110[1] = (this->_densityMapCorCoords[posIter])[1] * this->_ySamplingRate;
231  c110[2] = (this->_densityMapCorCoords[posIter])[2] * this->_zSamplingRate;
232  c110[3] = this->_densityMapCor[posIter];
233 
234  posIter = zTop + (this->_maxMapW + 1) * ( yTop + (this->_maxMapV + 1) * xTop );
235  if ( posIter > ( (this->_maxMapU+1) * (this->_maxMapV+1) * (this->_maxMapW+1) ) ) { this->_shellMappedData[sh][static_cast<unsigned int> ( phIt * this->_thetaAngle + thIt)] = 0.0; continue; }
236  c111[0] = (this->_densityMapCorCoords[posIter])[0] * this->_xSamplingRate;
237  c111[1] = (this->_densityMapCorCoords[posIter])[1] * this->_ySamplingRate;
238  c111[2] = (this->_densityMapCorCoords[posIter])[2] * this->_zSamplingRate;
239  c111[3] = this->_densityMapCor[posIter];
240 
241  //============================ Interpolate along x
242  xd = (x - ( ( xBottom - static_cast<int> ( ((this->_maxMapU+1)/2) ) ) * this->_xSamplingRate ) ) / this->_xSamplingRate;
243  c00[0] = (this->_xSamplingRate * xd) + c000[0]; c00[1] = c000[1]; c00[2] = c000[2]; c00[3] = ( c000[3] * ( 1.0 - xd ) ) + ( c100[3] * xd );
244  c01[0] = (this->_xSamplingRate * xd) + c001[0]; c01[1] = c001[1]; c01[2] = c001[2]; c01[3] = ( c001[3] * ( 1.0 - xd ) ) + ( c101[3] * xd );
245  c10[0] = (this->_xSamplingRate * xd) + c010[0]; c10[1] = c010[1]; c10[2] = c010[2]; c10[3] = ( c010[3] * ( 1.0 - xd ) ) + ( c110[3] * xd );
246  c11[0] = (this->_xSamplingRate * xd) + c011[0]; c11[1] = c011[1]; c11[2] = c011[2]; c11[3] = ( c011[3] * ( 1.0 - xd ) ) + ( c111[3] * xd );
247 
248  //============================ Interpolate along y
249  yd = (y - ( ( yBottom - static_cast<int> ( ((this->_maxMapV+1)/2) ) ) * this->_ySamplingRate ) ) / this->_ySamplingRate;
250  c0[0] = c00[0]; c0[1] = (this->_ySamplingRate * yd) + c00[1]; c0[2] = c00[2]; c0[3] = ( c00[3] * ( 1.0 - yd ) ) + ( c10[3] * yd );
251  c1[0] = c01[0]; c1[1] = (this->_ySamplingRate * yd) + c01[1]; c1[2] = c01[2]; c1[3] = ( c01[3] * ( 1.0 - yd ) ) + ( c11[3] * yd );
252 
253  //============================ Interpolate along z
254  zd = (z - ( ( zBottom - static_cast<int> ( ((this->_maxMapW+1)/2) ) ) * this->_zSamplingRate ) ) / this->_zSamplingRate;
255 
256  //============================ Save results
257  this->_shellMappedData[sh][static_cast<unsigned int> ( phIt * this->_thetaAngle + thIt)] = ( c0[3] * ( 1.0 - zd ) ) + ( c1[3] * zd );
258  }
259  }
260  }
261 
262  //======================================== Free unnecessary memory
263  delete[] this->_densityMapCorCoords;
264  this->_densityMapCorCoords = nullptr;
265 
266  if ( !keepInMemory )
267  {
268  delete[] this->_densityMapCor;
269  this->_densityMapCor = nullptr;
270  }
271 
272  //======================================== Done
273  this->_sphereMapped = true;
274 
275 }
276 
289 {
290  //======================================== Sanity check
291  if ( !this->_sphereMapped )
292  {
293  std::cerr << "!!! ProSHADE ERROR !!! Error in file " << this->_inputFileName << " !!! Attempted to obtain spherical harmonics before mapping the map onto spheres. Call the mapPhaselessToSphere function BEFORE the getSphericalHarmonicsCoeffs function." << std::endl;
294 
295  if ( settings->htmlReport )
296  {
297  std::stringstream hlpSS;
298  hlpSS << "<font color=\"red\">" << "Attempted to obtain spherical harmonics before mapping the map onto spheres. This looks like an internal bug, please report this case." << "</font>";
299  rvapi_set_text ( hlpSS.str().c_str(),
300  "ProgressSection",
301  settings->htmlReportLineProgress,
302  1,
303  1,
304  1 );
305  settings->htmlReportLineProgress += 1;
306  rvapi_flush ( );
307  }
308 
309  exit ( -1 );
310  }
311 
312  //======================================== Initialise internal values
313  this->_bandwidthLimit = bandwidth;
314  this->_oneDimmension = this->_bandwidthLimit * 2;
315 
316  //======================================== Allocate Memory
317  double *inputReal = new double [static_cast<unsigned int> ( this->_oneDimmension * this->_oneDimmension )];
318  double *inputZeroes = new double [static_cast<unsigned int> ( this->_oneDimmension * this->_oneDimmension )];
319  this->_realSHCoeffs = new double*[static_cast<unsigned int> ( this->_noShellsWithData )];
320  this->_imagSHCoeffs = new double*[static_cast<unsigned int> ( this->_noShellsWithData )];
321  this->_sphericalHarmonicsWeights = new double [static_cast<unsigned int> ( this->_bandwidthLimit * 4 )];
322  this->_semiNaiveTableSpace = new double [static_cast<unsigned int> ( Reduced_Naive_TableSize ( this->_bandwidthLimit, this->_bandwidthLimit ) +
323  Reduced_SpharmonicTableSize ( this->_bandwidthLimit, this->_bandwidthLimit ) )];
324 
325  this->_shWorkspace = (fftw_complex *) fftw_malloc ( sizeof(fftw_complex) * ( ( 8 * this->_bandwidthLimit * this->_bandwidthLimit ) +
326  ( 10 * this->_bandwidthLimit ) ) );
327 
328  for ( unsigned int shIt = 0; shIt < this->_noShellsWithData; shIt++ )
329  {
330  this->_realSHCoeffs[shIt] = new double [static_cast<unsigned int> ( this->_oneDimmension * this->_oneDimmension )];
331  this->_imagSHCoeffs[shIt] = new double [static_cast<unsigned int> ( this->_oneDimmension * this->_oneDimmension )];
332  }
333 
334  //======================================== Check Memory Allocation
335  if ( this->_shWorkspace == nullptr )
336  {
337  std::cerr << "!!! ProSHADE ERROR !!! Error in file " << this->_inputFileName << " !!! malloc (Memory Allocation) failed for _shWorkspace object in function getSphericalHarmonicsCoeffs." << std::endl;
338 
339  if ( settings->htmlReport )
340  {
341  std::stringstream hlpSS;
342  hlpSS << "<font color=\"red\">" << "Cannot allocate memory for spherical harmonics decomposition. Could you have run out of memory?" << "</font>";
343  rvapi_set_text ( hlpSS.str().c_str(),
344  "ProgressSection",
345  settings->htmlReportLineProgress,
346  1,
347  1,
348  1 );
349  settings->htmlReportLineProgress += 1;
350  rvapi_flush ( );
351  }
352 
353  exit( -1 ) ;
354  }
355 
356  //======================================== Generate Seminaive and naive table
357  this->_semiNaiveTable = SemiNaive_Naive_Pml_Table ( this->_bandwidthLimit,
358  this->_bandwidthLimit,
359  this->_semiNaiveTableSpace,
360  reinterpret_cast<double*> ( this->_shWorkspace ) );
361 
362  //======================================== Make weights for spherical transform
363  makeweights ( this->_bandwidthLimit,
364  this->_sphericalHarmonicsWeights );
365 
366  //======================================== Declare pointers to within workspace
367  double *rres;
368  double *ires;
369  double *fltres;
370  double *scratchpad;
371  double *rdataptr;
372  double *idataptr;
373 
374  //======================================== Initialize pointers to within workspace
375  rres = reinterpret_cast<double*> ( this->_shWorkspace );
376  ires = rres + ( this->_oneDimmension * this->_oneDimmension );
377  fltres = ires + ( this->_oneDimmension * this->_oneDimmension );
378  scratchpad = fltres + ( this->_bandwidthLimit );
379 
380  //======================================== Initialize fft plan along phi angles
381  fftw_plan fftPlan;
382  fftw_plan dctPlan;
383 
384  fftw_iodim dims[1];
385  fftw_iodim howmany_dims[1];
386 
387  int rank = 1;
388  int howmany_rank = 1;
389  double normCoeff = ( 1.0 / ( static_cast<double> ( this->_oneDimmension ) ) ) * sqrt( 2.0 * M_PI );
390  double powerOne = 1.0;
391  unsigned int hlp1 = 0;
392  unsigned int hlp2 = 0;
393 
394  dims[0].n = this->_oneDimmension;
395  dims[0].is = 1;
396  dims[0].os = this->_oneDimmension;
397 
398  howmany_dims[0].n = this->_oneDimmension;
399  howmany_dims[0].is = this->_oneDimmension;
400  howmany_dims[0].os = 1;
401 
402  //======================================== Plan fft transform
403  fftPlan = fftw_plan_guru_split_dft ( rank,
404  dims,
405  howmany_rank,
406  howmany_dims,
407  inputReal,
408  inputZeroes,
409  rres,
410  ires,
411  FFTW_ESTIMATE );
412 
413  //======================================== Initialize dct plan for SHT
414  dctPlan = fftw_plan_r2r_1d ( this->_oneDimmension,
415  scratchpad,
416  scratchpad + this->_oneDimmension,
417  FFTW_REDFT10,
418  FFTW_ESTIMATE ) ;
419 
420  //======================================== For each shell
421  for ( unsigned int shIter = 0; shIter < static_cast<unsigned int> ( this->_noShellsWithData ); shIter++ )
422  {
423  //==================================== Load pixel hits to signal array
424  for ( unsigned int iter = 0; iter < static_cast<unsigned int> ( this->_oneDimmension * this->_oneDimmension ); iter++ )
425  {
426  inputReal[iter] = this->_shellMappedData[shIter][iter];
427  inputZeroes[iter] = 0.0;
428  }
429 
430  //==================================== Execute fft plan along phi
431  fftw_execute_split_dft ( fftPlan, inputReal, inputZeroes, rres, ires ) ;
432 
433  //==================================== Normalize
434  for ( unsigned int iter = 0; iter < ( this->_oneDimmension * this->_oneDimmension ); iter++ )
435  {
436  rres[iter] *= normCoeff;
437  ires[iter] *= normCoeff;
438  }
439 
440  //==================================== Calculate the coefficients for each band
441  rdataptr = this->_realSHCoeffs[shIter];
442  idataptr = this->_imagSHCoeffs[shIter];
443  for ( unsigned int bandIter = 0; bandIter < this->_bandwidthLimit; bandIter++ )
444  {
445  //================================ *real* part calculation
446  SemiNaiveReduced ( rres + ( bandIter * this->_oneDimmension ),
447  this->_bandwidthLimit,
448  bandIter,
449  fltres,
450  scratchpad,
451  this->_semiNaiveTable[bandIter],
452  this->_sphericalHarmonicsWeights,
453  &dctPlan);
454 
455  //================================ Save the *real* results
456  memcpy ( rdataptr, fltres, sizeof(double) * ( this->_bandwidthLimit - bandIter ) );
457  rdataptr += this->_bandwidthLimit - bandIter;
458 
459  //================================ *imaginary* part calculation
460  SemiNaiveReduced ( ires + ( bandIter * this->_oneDimmension ),
461  this->_bandwidthLimit,
462  bandIter,
463  fltres,
464  scratchpad,
465  this->_semiNaiveTable[bandIter],
466  this->_sphericalHarmonicsWeights,
467  &dctPlan );
468 
469  //================================ Save the *imaginary* results
470  memcpy ( idataptr, fltres, sizeof(double) * ( this->_bandwidthLimit - bandIter ) );
471  idataptr += this->_bandwidthLimit - bandIter;
472  }
473 
474  //==================================== Since the data are real, there is no need to calculate negative order data, as there is the equivalence
475  powerOne = 1.0;
476  hlp1 = 0;
477  hlp2 = 0;
478  for( unsigned int order = 1; order < this->_bandwidthLimit; order++)
479  {
480  powerOne *= -1.0;
481  for( unsigned int bandIter = order; bandIter < this->_bandwidthLimit; bandIter++){
482 
483  hlp1 = seanindex ( order, bandIter, this->_bandwidthLimit );
484  hlp2 = seanindex ( -order, bandIter, this->_bandwidthLimit );
485 
486  this->_realSHCoeffs[shIter][hlp2] = powerOne * this->_realSHCoeffs[shIter][hlp1];
487  this->_imagSHCoeffs[shIter][hlp2] = -powerOne * this->_imagSHCoeffs[shIter][hlp1];
488  }
489  }
490  }
491 
492  delete[] inputReal;
493  delete[] inputZeroes;
494 
495  delete[] this->_semiNaiveTable;
496  delete[] this->_semiNaiveTableSpace;
497  delete[] this->_sphericalHarmonicsWeights;
498 
499  fftw_free ( this->_shWorkspace );
500 
501  this->_semiNaiveTable = nullptr;
502  this->_semiNaiveTableSpace = nullptr;
503  this->_sphericalHarmonicsWeights = nullptr;
504  this->_shWorkspace = nullptr;
505 
506  fftw_destroy_plan ( dctPlan );
507  fftw_destroy_plan ( fftPlan );
508 
509  //======================================== Done
510  this->_sphericalCoefficientsComputed = true;
511 
512 }
void mapPhaselessToSphere(ProSHADE::ProSHADE_settings *settings, double theta, double phi, double shellSz, unsigned int manualShells=0, bool keepInMemory=false, bool rotDefaults=false)
This function assumes the data have been processed and maps them onto a set of concentric spheres wit...
bool htmlReport
Should HTML report for the run be created?
Definition: ProSHADE.h:185
The main header file containing all declarations the user of the library needs.
void getSphericalHarmonicsCoeffs(unsigned int bandwidth, ProSHADE::ProSHADE_settings *settings)
This function takes the sphere mapped data and computes spherical harmoncis decomposition for each sh...
The main header file containing all declarations for the innter workings of the library.
int htmlReportLineProgress
Iterator for current HTML line in the progress bar.
Definition: ProSHADE.h:187
This class stores all the settings and is passed to the executive classes instead of multitude of par...
Definition: ProSHADE.h:74