ProSHADE  0.7.6.2 (DEC 2021)
Protein Shape Detection
ProSHADE_overlay.cpp
Go to the documentation of this file.
1 
22 //==================================================== ProSHADE
23 #include "ProSHADE_overlay.hpp"
24 
36 {
37  //================================================ Report progress
38  ProSHADE_internal_messages::printProgressMessage ( settings->verbose, 1, "Starting rotation function computation.", settings->messageShift );
39 
40  //================================================ Compute un-weighted E matrices and their weights
41  ProSHADE_internal_distances::computeEMatrices ( obj2, this, settings );
42 
43  //================================================ Normalise E matrices by the magnitudes
44  ProSHADE_internal_distances::normaliseEMatrices ( obj2, this, settings );
45 
46  //================================================ Generate SO(3) coefficients
48 
49  //================================================ Compute the inverse SO(3) Fourier Transform (SOFT) on the newly computed coefficients
51 
52  //================================================ Report completion
53  ProSHADE_internal_messages::printProgressMessage ( settings->verbose, 2, "Rotation function obtained.", settings->messageShift );
54 
55  //================================================ Done
56  return ;
57 
58 }
59 
74 void ProSHADE_internal_overlay::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 )
75 {
76  //================================================ Read in the structures
77  staticStructure->readInStructure ( settings->inputFiles.at(0), 0, settings );
78  movingStructure->readInStructure ( settings->inputFiles.at(1), 1, settings );
79 
80  //================================================ Internal data processing (COM, norm, mask, extra space)
81  staticStructure->processInternalMap ( settings );
82  movingStructure->processInternalMap ( settings );
83 
84  //================================================ Map to sphere
85  staticStructure->mapToSpheres ( settings );
86  movingStructure->mapToSpheres ( settings );
87 
88  //================================================ Get spherical harmonics
89  staticStructure->computeSphericalHarmonics ( settings );
90  movingStructure->computeSphericalHarmonics ( settings );
91 
92  //================================================ Get the rotation function of the pair
93  movingStructure->getOverlayRotationFunction ( settings, staticStructure );
94 
95  //================================================ Get inverse SO(3) map top peak Euler angle values
97  std::min ( staticStructure->getMaxBand(), movingStructure->getMaxBand() ) * 2,
98  eulA, eulB, eulG, settings );
99 
100  //================================================ Done
101  return ;
102 
103 }
104 
123 void ProSHADE_internal_overlay::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 )
124 {
125  //================================================ Read in the structures
126  staticStructure->readInStructure ( settings->inputFiles.at(0), 0, settings );
127  movingStructure->readInStructure ( settings->inputFiles.at(1), 1, settings );
128 
129  //================================================ Determine spherical harmonics variables from the first structure (otherwise, they would be determined from the second structure)
130  settings->determineAllSHValues ( staticStructure->xDimIndices, staticStructure->yDimIndices,
131  staticStructure->xDimSize, staticStructure->yDimSize, staticStructure->zDimSize );
132 
133  //================================================ Internal data processing (COM, norm, mask, extra space)
134  staticStructure->processInternalMap ( settings );
135  movingStructure->processInternalMap ( settings );
136 
137  //================================================ Rotate map
138  std::vector< proshade_double > rotCen = movingStructure->rotateMapRealSpaceInPlace ( eulA, eulB, eulG );
139 
140  //================================================ Zero padding for smaller structure
141  staticStructure->zeroPaddToDims ( std::max ( staticStructure->getXDim(), movingStructure->getXDim() ),
142  std::max ( staticStructure->getYDim(), movingStructure->getYDim() ),
143  std::max ( staticStructure->getZDim(), movingStructure->getZDim() ) );
144  movingStructure->zeroPaddToDims ( std::max ( staticStructure->getXDim(), movingStructure->getXDim() ),
145  std::max ( staticStructure->getYDim(), movingStructure->getYDim() ),
146  std::max ( staticStructure->getZDim(), movingStructure->getZDim() ) );
147 
148  //================================================ Report progress
149  ProSHADE_internal_messages::printProgressMessage ( settings->verbose, 1, "Starting translation function computation." );
150 
151  //================================================ Compute the translation function map
152  movingStructure->computeTranslationMap ( staticStructure );
153 
154  //================================================ Find highest peak in translation function
155  proshade_double mapPeak = 0.0;
156  proshade_unsign xDimS = staticStructure->getXDim();
157  proshade_unsign yDimS = staticStructure->getYDim();
158  proshade_unsign zDimS = staticStructure->getZDim();
159  ProSHADE_internal_maths::findHighestValueInMap ( movingStructure->getTranslationFnPointer(), xDimS, yDimS, zDimS, trsX, trsY, trsZ, &mapPeak );
160 
161  //================================================ Dont translate over half
162  if ( *trsX > ( static_cast< proshade_double > ( xDimS ) / 2.0 ) ) { *trsX = *trsX - static_cast< proshade_double > ( xDimS ); }
163  if ( *trsY > ( static_cast< proshade_double > ( yDimS ) / 2.0 ) ) { *trsY = *trsY - static_cast< proshade_double > ( yDimS ); }
164  if ( *trsZ > ( static_cast< proshade_double > ( zDimS ) / 2.0 ) ) { *trsZ = *trsZ - static_cast< proshade_double > ( zDimS ); }
165 
166  //================================================ Convert the map positions onto translation in Angstroms and save the results
167  computeTranslationsFromPeak ( staticStructure, movingStructure, trsX, trsY, trsZ );
168 
169  movingStructure->originalPdbTransX = *trsX;
170  movingStructure->originalPdbTransY = *trsY;
171  movingStructure->originalPdbTransZ = *trsZ;
172 
173  proshade_single xMov = static_cast< proshade_single > ( -(*trsX) );
174  proshade_single yMov = static_cast< proshade_single > ( -(*trsY) );
175  proshade_single zMov = static_cast< proshade_single > ( -(*trsZ) );
176 
177  //================================================ Report progress
178  std::stringstream hlpSS;
179  hlpSS << "Optimal map translation distances are " << *trsX << " ; " << *trsY << " ; " << *trsZ << " Angstroms with peak height " << mapPeak / ( static_cast< proshade_double > ( xDimS ) * static_cast< proshade_double > ( yDimS ) * static_cast< proshade_double > ( zDimS ) );
180 
181  //================================================ Save original from variables for PDB output
182  movingStructure->mapMovFromsChangeX = static_cast< proshade_double > ( movingStructure->xFrom );
183  movingStructure->mapMovFromsChangeY = static_cast< proshade_double > ( movingStructure->yFrom );
184  movingStructure->mapMovFromsChangeZ = static_cast< proshade_double > ( movingStructure->zFrom );
185 
186  //================================================ Move the map
187  ProSHADE_internal_mapManip::moveMapByIndices ( &xMov, &yMov, &zMov, movingStructure->getXDimSize(), movingStructure->getYDimSize(), movingStructure->getZDimSize(),
188  movingStructure->getXFromPtr(), movingStructure->getXToPtr(),
189  movingStructure->getYFromPtr(), movingStructure->getYToPtr(),
190  movingStructure->getZFromPtr(), movingStructure->getZToPtr(),
191  movingStructure->getXAxisOrigin(), movingStructure->getYAxisOrigin(), movingStructure->getZAxisOrigin() );
192 
193  ProSHADE_internal_mapManip::moveMapByFourier ( movingStructure->getInternalMap(), xMov, yMov, zMov,
194  movingStructure->getXDimSize(), movingStructure->getYDimSize(), movingStructure->getZDimSize(),
195  static_cast< proshade_signed > ( movingStructure->getXDim() ), static_cast< proshade_signed > ( movingStructure->getYDim() ),
196  static_cast< proshade_signed > ( movingStructure->getZDim() ) );
197 
198  //================================================ Report progress
199  ProSHADE_internal_messages::printProgressMessage ( settings->verbose, 3, hlpSS.str(), settings->messageShift );
200  ProSHADE_internal_messages::printProgressMessage ( settings->verbose, 2, "Translation function computation complete.", settings->messageShift );
201 
202  //================================================ Keep only the change in from and to variables
203  movingStructure->mapMovFromsChangeX = static_cast< proshade_double > ( movingStructure->xFrom ) - movingStructure->mapMovFromsChangeX;
204  movingStructure->mapMovFromsChangeY = static_cast< proshade_double > ( movingStructure->yFrom ) - movingStructure->mapMovFromsChangeY;
205  movingStructure->mapMovFromsChangeZ = static_cast< proshade_double > ( movingStructure->zFrom ) - movingStructure->mapMovFromsChangeZ;
206 
207  //================================================ Done
208  return ;
209 
210 }
211 
220 void ProSHADE_internal_overlay::computeTranslationsFromPeak ( ProSHADE_internal_data::ProSHADE_data* staticStructure, ProSHADE_internal_data::ProSHADE_data* movingStructure, proshade_double *trsX, proshade_double *trsY, proshade_double *trsZ )
221 {
222  //================================================ Compute map movement
223  proshade_single xSamplRate = ( static_cast< proshade_single > ( staticStructure->getXDimSize() ) / static_cast< proshade_single > ( staticStructure->getXDim() ) );
224  proshade_single ySamplRate = ( static_cast< proshade_single > ( staticStructure->getYDimSize() ) / static_cast< proshade_single > ( staticStructure->getYDim() ) );
225  proshade_single zSamplRate = ( static_cast< proshade_single > ( staticStructure->getZDimSize() ) / static_cast< proshade_single > ( staticStructure->getZDim() ) );
226 
227  proshade_single xCor = static_cast< proshade_single > ( static_cast< proshade_single > ( staticStructure->xFrom - movingStructure->xFrom ) * xSamplRate );
228  proshade_single yCor = static_cast< proshade_single > ( static_cast< proshade_single > ( staticStructure->yFrom - movingStructure->yFrom ) * ySamplRate );
229  proshade_single zCor = static_cast< proshade_single > ( static_cast< proshade_single > ( staticStructure->zFrom - movingStructure->zFrom ) * zSamplRate );
230 
231  proshade_single xMov = static_cast< proshade_single > ( 0.0 - static_cast< proshade_double > ( xCor ) - ( *trsX * static_cast< proshade_double > ( xSamplRate ) ) );
232  proshade_single yMov = static_cast< proshade_single > ( 0.0 - static_cast< proshade_double > ( yCor ) - ( *trsY * static_cast< proshade_double > ( ySamplRate ) ) );
233  proshade_single zMov = static_cast< proshade_single > ( 0.0 - static_cast< proshade_double > ( zCor ) - ( *trsZ * static_cast< proshade_double > ( zSamplRate ) ) );
234 
235  //================================================ Save translation vector back
236  *trsX = static_cast< proshade_double > ( -xMov );
237  *trsY = static_cast< proshade_double > ( -yMov );
238  *trsZ = static_cast< proshade_double > ( -zMov );
239 
240  //================================================ Done
241  return ;
242 
243 }
244 
251 {
252  //================================================ Initialise local variables
253  std::vector< proshade_double > ret;
254  proshade_double mapPeak = 0.0;
255  proshade_unsign xDimS = staticStructure->getXDim();
256  proshade_unsign yDimS = staticStructure->getYDim();
257  proshade_unsign zDimS = staticStructure->getZDim();
258  proshade_double trsX = 0.0, trsY = 0.0, trsZ = 0.0;
259  ProSHADE_internal_maths::findHighestValueInMap ( this->getTranslationFnPointer(),
260  xDimS,
261  yDimS,
262  zDimS,
263  &trsX,
264  &trsY,
265  &trsZ,
266  &mapPeak );
267 
268  //================================================ Dont translate over half
269  if ( trsX > ( static_cast< proshade_double > ( xDimS ) / 2.0 ) ) { trsX = trsX - static_cast< proshade_double > ( xDimS ); }
270  if ( trsY > ( static_cast< proshade_double > ( yDimS ) / 2.0 ) ) { trsY = trsY - static_cast< proshade_double > ( yDimS ); }
271  if ( trsZ > ( static_cast< proshade_double > ( zDimS ) / 2.0 ) ) { trsZ = trsZ - static_cast< proshade_double > ( zDimS ); }
272 
273  //================================================ Convert the map positions onto translation in Angstroms and save the results
274  ProSHADE_internal_overlay::computeTranslationsFromPeak ( staticStructure, this, &trsX, &trsY, &trsZ );
275 
276  this->originalPdbTransX = trsX;
277  this->originalPdbTransY = trsY;
278  this->originalPdbTransZ = trsZ;
279 
280  proshade_single xMov = static_cast< proshade_single > ( -trsX );
281  proshade_single yMov = static_cast< proshade_single > ( -trsY );
282  proshade_single zMov = static_cast< proshade_single > ( -trsZ );
283 
284  //================================================ Save results as vector
285  ProSHADE_internal_misc::addToDoubleVector ( &ret, static_cast<proshade_double> ( -xMov ) );
286  ProSHADE_internal_misc::addToDoubleVector ( &ret, static_cast<proshade_double> ( -yMov ) );
287  ProSHADE_internal_misc::addToDoubleVector ( &ret, static_cast<proshade_double> ( -zMov ) );
288 
289  //================================================ Save original from variables for PDB output
290  this->mapMovFromsChangeX = static_cast<proshade_double> ( this->xFrom );
291  this->mapMovFromsChangeY = static_cast<proshade_double> ( this->yFrom );
292  this->mapMovFromsChangeZ = static_cast<proshade_double> ( this->zFrom );
293 
294  //================================================ Move the map
295  ProSHADE_internal_mapManip::moveMapByIndices ( &xMov, &yMov, &zMov, this->getXDimSize(), this->getYDimSize(), this->getZDimSize(),
296  this->getXFromPtr(), this->getXToPtr(),
297  this->getYFromPtr(), this->getYToPtr(),
298  this->getZFromPtr(), this->getZToPtr(),
299  this->getXAxisOrigin(), this->getYAxisOrigin(), this->getZAxisOrigin() );
300 
301  ProSHADE_internal_mapManip::moveMapByFourier ( this->getInternalMap(), xMov, yMov, zMov,
302  this->getXDimSize(), this->getYDimSize(), this->getZDimSize(),
303  static_cast< proshade_signed > ( this->getXDim() ), static_cast< proshade_signed > ( this->getYDim() ), static_cast< proshade_signed > ( this->getZDim() ) );
304 
305  //================================================ Keep only the change in from and to variables
306  this->mapMovFromsChangeX = static_cast<proshade_double> ( this->xFrom ) - this->mapMovFromsChangeX;
307  this->mapMovFromsChangeY = static_cast<proshade_double> ( this->yFrom ) - this->mapMovFromsChangeY;
308  this->mapMovFromsChangeZ = static_cast<proshade_double> ( this->zFrom ) - this->mapMovFromsChangeZ;
309 
310  //================================================ Done
311  return ( ret );
312 
313 }
314 
325 {
326  //================================================ Do this using Fourier!
327  fftw_complex *tmpIn1 = nullptr, *tmpOut1 = nullptr, *tmpIn2 = nullptr, *tmpOut2 = nullptr, *resOut = nullptr;
328  fftw_plan forwardFourierObj1, forwardFourierObj2, inverseFourierCombo;
329  proshade_unsign dimMult = staticStructure->getXDim() * staticStructure->getYDim() * staticStructure->getZDim();
330  ProSHADE_internal_overlay::allocateTranslationFunctionMemory ( tmpIn1, tmpOut1, tmpIn2, tmpOut2, this->translationMap, resOut, forwardFourierObj1, forwardFourierObj2, inverseFourierCombo, staticStructure->getXDim(), staticStructure->getYDim(), staticStructure->getZDim() );
331 
332  //================================================ Fill in input data
333  for ( proshade_unsign iter = 0; iter < dimMult; iter++ ) { tmpIn1[iter][0] = staticStructure->getMapValue ( iter ); tmpIn1[iter][1] = 0.0; }
334  for ( proshade_unsign iter = 0; iter < dimMult; iter++ ) { tmpIn2[iter][0] = this->getMapValue ( iter ); tmpIn2[iter][1] = 0.0; }
335 
336  //================================================ Calculate Fourier
337  fftw_execute ( forwardFourierObj1 );
338  fftw_execute ( forwardFourierObj2 );
339 
340  //================================================ Combine Fourier coeffs and invert
341  ProSHADE_internal_maths::combineFourierForTranslation ( tmpOut1, tmpOut2, resOut, staticStructure->getXDim(), staticStructure->getYDim(), staticStructure->getZDim() );
342  fftw_execute ( inverseFourierCombo );
343 
344  //================================================ Free memory
345  ProSHADE_internal_overlay::freeTranslationFunctionMemory ( tmpIn1, tmpOut1, tmpIn2, tmpOut2, resOut, forwardFourierObj1, forwardFourierObj2, inverseFourierCombo );
346 
347  //================================================ Done
348  return ;
349 
350 }
351 
367 void ProSHADE_internal_overlay::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 )
368 {
369  //================================================ Allocate memory
370  tmpIn1 = new fftw_complex[xD * yD * zD];
371  tmpOut1 = new fftw_complex[xD * yD * zD];
372  tmpIn2 = new fftw_complex[xD * yD * zD];
373  tmpOut2 = new fftw_complex[xD * yD * zD];
374  resIn = new fftw_complex[xD * yD * zD];
375  resOut = new fftw_complex[xD * yD * zD];
376 
377  //================================================ Check memory allocation
378  ProSHADE_internal_misc::checkMemoryAllocation ( tmpIn1 , __FILE__, __LINE__, __func__ );
379  ProSHADE_internal_misc::checkMemoryAllocation ( tmpOut1, __FILE__, __LINE__, __func__ );
380  ProSHADE_internal_misc::checkMemoryAllocation ( tmpIn2, __FILE__, __LINE__, __func__ );
381  ProSHADE_internal_misc::checkMemoryAllocation ( tmpOut2, __FILE__, __LINE__, __func__ );
382  ProSHADE_internal_misc::checkMemoryAllocation ( resIn, __FILE__, __LINE__, __func__ );
383  ProSHADE_internal_misc::checkMemoryAllocation ( resOut, __FILE__, __LINE__, __func__ );
384 
385  //================================================ Get Fourier transforms of the maps
386  forwardFourierObj1 = fftw_plan_dft_3d ( static_cast< int > ( xD ), static_cast< int > ( yD ), static_cast< int > ( zD ), tmpIn1, tmpOut1, FFTW_FORWARD , FFTW_ESTIMATE );
387  forwardFourierObj2 = fftw_plan_dft_3d ( static_cast< int > ( xD ), static_cast< int > ( yD ), static_cast< int > ( zD ), tmpIn2, tmpOut2, FFTW_FORWARD , FFTW_ESTIMATE );
388  inverseFourierCombo = fftw_plan_dft_3d ( static_cast< int > ( xD ), static_cast< int > ( yD ), static_cast< int > ( zD ), resOut, resIn , FFTW_BACKWARD, FFTW_ESTIMATE );
389 
390  //================================================ Done
391  return ;
392 
393 }
394 
406 void ProSHADE_internal_overlay::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 )
407 {
408  //================================================ Release memory
409  fftw_destroy_plan ( forwardFourierObj1 );
410  fftw_destroy_plan ( forwardFourierObj2 );
411  fftw_destroy_plan ( inverseFourierCombo );
412  delete[] tmpIn1;
413  delete[] tmpIn2;
414  delete[] tmpOut1;
415  delete[] tmpOut2;
416  delete[] resOut;
417 
418  //======================================== Done
419  return ;
420 
421 }
422 
433 void ProSHADE_internal_data::ProSHADE_data::zeroPaddToDims ( proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim )
434 {
435  //================================================ Sanity check
436  if ( ( this->xDimIndices > xDim ) || ( this->yDimIndices > yDim ) || ( this->zDimIndices > zDim ) )
437  {
438  throw ProSHADE_exception ( "Cannot zero-pad in negative direction.", "EO00034", __FILE__, __LINE__, __func__, "The requested padded size of a structure is smaller than\n : the current size. If the user sees this error, there is\n : likely a considerable bug. Please report this error." );
439  }
440 
441  //================================================ If done, do nothing
442  if ( ( this->xDimIndices == xDim ) && ( this->yDimIndices == yDim ) && ( this->zDimIndices == zDim ) ) { return ; }
443 
444  //================================================ Find out how many zeroes need to be added before and after
445  proshade_unsign addXPre, addYPre, addZPre, addXPost, addYPost, addZPost;
446  ProSHADE_internal_overlay::computeBeforeAfterZeroCounts ( &addXPre, &addYPre, &addZPre, &addXPost, &addYPost, &addZPost, xDim, yDim, zDim, this->xDimIndices, this->yDimIndices, this->zDimIndices );
447 
448  //================================================ Create a new map
449  proshade_double* newMap = new proshade_double [xDim * yDim * zDim];
450 
451  //================================================ Do the hard work
452  ProSHADE_internal_overlay::paddMapWithZeroes ( this->internalMap, newMap, xDim, yDim, zDim, this->xDimIndices, this->yDimIndices, this->zDimIndices, addXPre, addYPre, addZPre );
453 
454  //================================================ Create a new internal map and copy
455  delete[] this->internalMap;
456  this->internalMap = new proshade_double [xDim * yDim * zDim];
457  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( xDim * yDim * zDim ); iter++ ) { this->internalMap[iter] = newMap[iter]; }
458 
459  //================================================ Release memory
460  delete[] newMap;
461 
462  //================================================ Update map related variables
463  this->xDimSize = static_cast< proshade_single > ( xDim ) * ( this->xDimSize / static_cast< proshade_single > ( this->xDimIndices ) );
464  this->yDimSize = static_cast< proshade_single > ( yDim ) * ( this->yDimSize / static_cast< proshade_single > ( this->yDimIndices ) );
465  this->zDimSize = static_cast< proshade_single > ( zDim ) * ( this->zDimSize / static_cast< proshade_single > ( this->zDimIndices ) );
466  this->xDimIndices = xDim ; this->yDimIndices = yDim ; this->zDimIndices = zDim;
467  this->xGridIndices = xDim ; this->yGridIndices = yDim ; this->zGridIndices = zDim;
468  this->xFrom -= addXPre ; this->yFrom -= addYPre ; this->zFrom -= addZPre;
469  this->xTo += addXPost; this->yTo += addYPost; this->zTo += addZPost;
470  this->xAxisOrigin -= addXPre ; this->yAxisOrigin -= addYPre ; this->zAxisOrigin -= addZPre ;
471 
472  //================================================ Done
473  return ;
474 
475 }
476 
492 void ProSHADE_internal_overlay::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 )
493 {
494  //================================================ Compute
495  *addXPre = ( xDim - xDimIndices ) / 2;
496  *addYPre = ( yDim - yDimIndices ) / 2;
497  *addZPre = ( zDim - zDimIndices ) / 2;
498  *addXPost = *addXPre; if ( ( xDim - xDimIndices ) % 2 == 1 ) { *addXPost += 1; }
499  *addYPost = *addYPre; if ( ( yDim - yDimIndices ) % 2 == 1 ) { *addYPost += 1; }
500  *addZPost = *addZPre; if ( ( zDim - zDimIndices ) % 2 == 1 ) { *addZPost += 1; }
501 
502  //================================================ Done
503  return ;
504 
505 }
506 
521 void ProSHADE_internal_overlay::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 )
522 {
523  //================================================ Initialise local variables
524  proshade_unsign newMapIndex = 0;
525  proshade_unsign oldMapIndex = 0;
526 
527  //================================================ Copy and padd as appropriate
528  for ( proshade_unsign xIt = 0; xIt < xDim; xIt++ )
529  {
530  for ( proshade_unsign yIt = 0; yIt < yDim; yIt++ )
531  {
532  for ( proshade_unsign zIt = 0; zIt < zDim; zIt++ )
533  {
534  //==================================== Find map location
535  newMapIndex = zIt + zDim * ( yIt + yDim * xIt );
536 
537  //==================================== If less than needed, add zero
538  if ( xIt < addXPre ) { newMap[newMapIndex] = 0.0; continue; }
539  if ( yIt < addYPre ) { newMap[newMapIndex] = 0.0; continue; }
540  if ( zIt < addZPre ) { newMap[newMapIndex] = 0.0; continue; }
541 
542  //==================================== If more than needed, add zero
543  if ( xIt >= ( addXPre + xDimIndices ) ) { newMap[newMapIndex] = 0.0; continue; }
544  if ( yIt >= ( addYPre + yDimIndices ) ) { newMap[newMapIndex] = 0.0; continue; }
545  if ( zIt >= ( addZPre + zDimIndices ) ) { newMap[newMapIndex] = 0.0; continue; }
546 
547  //==================================== If not padding, copy original values
548  oldMapIndex = (zIt-addZPre) + zDimIndices * ( (yIt-addYPre) + yDimIndices * (xIt-addXPre) );
549  newMap[newMapIndex] = oldMap[oldMapIndex];
550  }
551  }
552  }
553 
554  //======================================== Done
555  return ;
556 
557 }
558 
571 void ProSHADE_internal_data::ProSHADE_data::rotateMapReciprocalSpace ( ProSHADE_settings* settings, proshade_double eulerAlpha, proshade_double eulerBeta, proshade_double eulerGamma )
572 {
573  //================================================ Set maximum comparison bandwidth to maximum object bandwidth
574  this->maxCompBand = this->spheres[this->noSpheres-1]->getLocalBandwidth();
575 
576  //================================================ Save map COM after processing but before rotation
577  this->findMapCOM ( );
578  this->mapCOMProcessChangeX += ( this->xCom - this->originalMapXCom );
579  this->mapCOMProcessChangeY += ( this->yCom - this->originalMapYCom );
580  this->mapCOMProcessChangeZ += ( this->zCom - this->originalMapZCom );
581 
582  //================================================ Compute the Wigner D matrices for the Euler angles
583  ProSHADE_internal_wigner::computeWignerMatricesForRotation ( settings, this, -eulerAlpha, eulerBeta, -eulerGamma );
584 
585  //================================================ Initialise rotated Spherical Harmonics memory
586  this->allocateRotatedSHMemory ( );
587 
588  //================================================ Multiply SH coeffs by Wigner
589  this->computeRotatedSH ( );
590 
591  //================================================ Inverse the SH coeffs to shells
592  this->invertSHCoefficients ( );
593 
594  //================================================ Find spherical cut-offs
595  std::vector<proshade_double> lonCO, latCO;
596  ProSHADE_internal_overlay::computeAngularThreshold ( &lonCO, &latCO, settings->maxBandwidth * 2 );
597 
598  //================================================ Allocate memory for the rotated map
599  proshade_double *densityMapRotated = new proshade_double [this->xDimIndices * this->yDimIndices * this->zDimIndices];
600  ProSHADE_internal_misc::checkMemoryAllocation ( densityMapRotated, __FILE__, __LINE__, __func__ );
601  for ( unsigned int iter = 0; iter < static_cast<unsigned int> ( this->xDimIndices * this->yDimIndices * this->zDimIndices ); iter++ ) { densityMapRotated[iter] = 0.0; }
602 
603  //================================================ Interpolate onto cartesian grid
604  this->interpolateMapFromSpheres ( densityMapRotated );
605 
606  //================================================ Copy map
607  for ( proshade_unsign iter = 0; iter < ( this->xDimIndices * this->yDimIndices * this->zDimIndices ); iter++ )
608  {
609  this->internalMap[iter] = densityMapRotated[iter];
610  }
611 
612  //================================================ Release rotated map (original is now rotated)
613  delete[] densityMapRotated;
614 
615  //================================================ Done
616  return ;
617 
618 }
619 
632 std::vector< proshade_double > ProSHADE_internal_data::ProSHADE_data::rotateMapRealSpace ( proshade_double axX, proshade_double axY, proshade_double axZ, proshade_double axAng, proshade_double*& map )
633 {
634  //================================================ Initialise local variables
635  bool withinBounds = true;
636  proshade_double c000, c001, c010, c011, c100, c101, c110, c111, c00, c01, c10, c11, c0, c1;
637  size_t arrPos = 0;
638  proshade_double xCOM, yCOM, zCOM;
639  std::vector< proshade_double > ret;
640 
641  //================================================ Store sampling rates
642  proshade_single xSampRate = this->xDimSize / static_cast< proshade_single > ( this->xTo - this->xFrom );
643  proshade_single ySampRate = this->yDimSize / static_cast< proshade_single > ( this->yTo - this->yFrom );
644  proshade_single zSampRate = this->zDimSize / static_cast< proshade_single > ( this->zTo - this->zFrom );
645 
646  //================================================ Compute map COM
647  ProSHADE_internal_mapManip::findMAPCOMValues ( this->internalMap, &xCOM, &yCOM, &zCOM, this->xDimSize, this->yDimSize, this->zDimSize, this->xFrom, this->xTo, this->yFrom, this->yTo, this->zFrom, this->zTo );
648 
649  //================================================ Allocate local variables
650  proshade_single *mins = new proshade_single[3];
651  proshade_single *maxs = new proshade_single[3];
652  proshade_single *rotMat = new proshade_single[9];
653  proshade_single *rotVec;
654  proshade_single *interpMins = new proshade_single[3];
655  proshade_single *interpMaxs = new proshade_single[3];
656  proshade_single *interpDiff = new proshade_single[3];
657  proshade_single *movs = new proshade_single[3];
658  map = new proshade_double[ this->xDimIndices * this->yDimIndices * this->zDimIndices ];
659 
660  //================================================ Check memory allocation
661  ProSHADE_internal_misc::checkMemoryAllocation ( mins, __FILE__, __LINE__, __func__ );
662  ProSHADE_internal_misc::checkMemoryAllocation ( maxs, __FILE__, __LINE__, __func__ );
663  ProSHADE_internal_misc::checkMemoryAllocation ( rotMat, __FILE__, __LINE__, __func__ );
664  ProSHADE_internal_misc::checkMemoryAllocation ( interpMins, __FILE__, __LINE__, __func__ );
665  ProSHADE_internal_misc::checkMemoryAllocation ( interpMaxs, __FILE__, __LINE__, __func__ );
666  ProSHADE_internal_misc::checkMemoryAllocation ( interpDiff, __FILE__, __LINE__, __func__ );
667  ProSHADE_internal_misc::checkMemoryAllocation ( movs, __FILE__, __LINE__, __func__ );
668  ProSHADE_internal_misc::checkMemoryAllocation ( map, __FILE__, __LINE__, __func__ );
669 
670  //================================================ Fill map with zeroes
671  for ( size_t iter = 0; iter < static_cast< size_t > ( this->xDimIndices * this->yDimIndices * this->zDimIndices ); iter++ ) { map[iter] = 0.0; }
672 
673  //================================================ Determine map max's and min's in terms of the hkl system
674  mins[0] = std::floor ( static_cast< proshade_single > ( this->xDimIndices ) / -2.0f );
675  mins[1] = std::floor ( static_cast< proshade_single > ( this->yDimIndices ) / -2.0f );
676  mins[2] = std::floor ( static_cast< proshade_single > ( this->zDimIndices ) / -2.0f );
677 
678  maxs[0] = -mins[0];
679  maxs[1] = -mins[1];
680  maxs[2] = -mins[2];
681 
682  if ( this->xDimIndices % 2 == 0 ) { maxs[0] -= 1.0f; }
683  if ( this->yDimIndices % 2 == 0 ) { maxs[1] -= 1.0f; }
684  if ( this->zDimIndices % 2 == 0 ) { maxs[2] -= 1.0f; }
685 
686  //================================================ Rotate about COM instead of map midpoint
687  movs[0] = 0.0; //( mins[0] - static_cast< proshade_single > ( this->xFrom ) );
688  movs[1] = 0.0; //( mins[1] - static_cast< proshade_single > ( this->yFrom ) );
689  movs[2] = 0.0; //( mins[2] - static_cast< proshade_single > ( this->zFrom ) );
690 
691  //================================================ Save rotation centre
692  ProSHADE_internal_misc::addToDoubleVector ( &ret, static_cast< proshade_double > ( ( -mins[0] * xSampRate ) + ( static_cast< proshade_single > ( this->xFrom ) * xSampRate ) ) );
693  ProSHADE_internal_misc::addToDoubleVector ( &ret, static_cast< proshade_double > ( ( -mins[1] * ySampRate ) + ( static_cast< proshade_single > ( this->yFrom ) * ySampRate ) ) );
694  ProSHADE_internal_misc::addToDoubleVector ( &ret, static_cast< proshade_double > ( ( -mins[2] * zSampRate ) + ( static_cast< proshade_single > ( this->zFrom ) * zSampRate ) ) );
695 
696  //================================================ Get rotation matrix from Euler angles
697  ProSHADE_internal_maths::getRotationMatrixFromAngleAxis ( rotMat, axX, axY, axZ, axAng );
698 
699  //================================================ For each point
700  for ( proshade_single xIt = mins[0]; xIt <= maxs[0]; xIt += 1.0f )
701  {
702  for ( proshade_single yIt = mins[1]; yIt <= maxs[1]; yIt += 1.0f )
703  {
704  for ( proshade_single zIt = mins[2]; zIt <= maxs[2]; zIt += 1.0f )
705  {
706  //==================================== Compute new point position
707  rotVec = ProSHADE_internal_maths::compute3x3MatrixVectorMultiplication ( rotMat, xIt - movs[0], yIt - movs[1], zIt - movs[2] );
708 
709  //==================================== Find surrounding grid points indices and check for boundaries
710  withinBounds = true;
711  for ( size_t posIt = 0; posIt < 3; posIt++ )
712  {
713  //================================ Determine surrounding points indices in this dimension
714  interpMins[posIt] = std::floor ( rotVec[posIt] );
715  interpMaxs[posIt] = interpMins[posIt] + 1.0f;
716 
717  //================================ Check for boundaries
718  if ( ( maxs[posIt] < interpMins[posIt] ) || ( interpMins[posIt] < mins[posIt] ) || ( maxs[posIt] < interpMaxs[posIt] ) || ( interpMaxs[posIt] < mins[posIt] ) )
719  {
720  withinBounds = false;
721  break;
722  }
723 
724  //================================ Compute the difference from position to min index along this axis
725  interpDiff[posIt] = rotVec[posIt] - interpMins[posIt];
726  }
727  if ( !withinBounds ) { continue; }
728 
729  //==================================== Make sure interpolation max's are within bounds
730  for ( size_t posIt = 0; posIt < 3; posIt++ )
731  {
732  interpMaxs[posIt] = std::min ( maxs[posIt], std::max ( mins[posIt], interpMaxs[posIt] ) );
733  }
734 
735  //==================================== Release memory
736  delete[] rotVec;
737 
738  //==================================== Find the surrounding points values from their indices
739  arrPos = static_cast< size_t > ( ( interpMins[2] - mins[2] ) + static_cast< proshade_single > ( this->zDimIndices ) * ( ( interpMins[1] - mins[1] ) + static_cast< proshade_single > ( this->yDimIndices ) * ( interpMins[0] - mins[0] ) ) );
740  c000 = this->internalMap[arrPos];
741 
742  arrPos = static_cast< size_t > ( ( interpMaxs[2] - mins[2] ) + static_cast< proshade_single > ( this->zDimIndices ) * ( ( interpMins[1] - mins[1] ) + static_cast< proshade_single > ( this->yDimIndices ) * ( interpMins[0] - mins[0] ) ) );
743  c001 = this->internalMap[arrPos];
744 
745  arrPos = static_cast< size_t > ( ( interpMins[2] - mins[2] ) + static_cast< proshade_single > ( this->zDimIndices ) * ( ( interpMaxs[1] - mins[1] ) + static_cast< proshade_single > ( this->yDimIndices ) * ( interpMins[0] - mins[0] ) ) );
746  c010 = this->internalMap[arrPos];
747 
748  arrPos = static_cast< size_t > ( ( interpMaxs[2] - mins[2] ) + static_cast< proshade_single > ( this->zDimIndices ) * ( ( interpMaxs[1] - mins[1] ) + static_cast< proshade_single > ( this->yDimIndices ) * ( interpMins[0] - mins[0] ) ) );
749  c011 = this->internalMap[arrPos];
750 
751  arrPos = static_cast< size_t > ( ( interpMins[2] - mins[2] ) + static_cast< proshade_single > ( this->zDimIndices ) * ( ( interpMins[1] - mins[1] ) + static_cast< proshade_single > ( this->yDimIndices ) * ( interpMaxs[0] - mins[0] ) ) );
752  c100 = this->internalMap[arrPos];
753 
754  arrPos = static_cast< size_t > ( ( interpMaxs[2] - mins[2] ) + static_cast< proshade_single > ( this->zDimIndices ) * ( ( interpMins[1] - mins[1] ) + static_cast< proshade_single > ( this->yDimIndices ) * ( interpMaxs[0] - mins[0] ) ) );
755  c101 = this->internalMap[arrPos];
756 
757  arrPos = static_cast< size_t > ( ( interpMins[2] - mins[2] ) + static_cast< proshade_single > ( this->zDimIndices ) * ( ( interpMaxs[1] - mins[1] ) + static_cast< proshade_single > ( this->yDimIndices ) * ( interpMaxs[0] - mins[0] ) ) );
758  c110 = this->internalMap[arrPos];
759 
760  arrPos = static_cast< size_t > ( ( interpMaxs[2] - mins[2] ) + static_cast< proshade_single > ( this->zDimIndices ) * ( ( interpMaxs[1] - mins[1] ) + static_cast< proshade_single > ( this->yDimIndices ) * ( interpMaxs[0] - mins[0] ) ) );
761  c111 = this->internalMap[arrPos];
762 
763  //==================================== Interpolate along x-axis
764  c00 = ( c000 * ( 1.0 - static_cast< proshade_double > ( interpDiff[0] ) ) ) + ( c100 * static_cast< proshade_double > ( interpDiff[0] ) );
765  c01 = ( c001 * ( 1.0 - static_cast< proshade_double > ( interpDiff[0] ) ) ) + ( c101 * static_cast< proshade_double > ( interpDiff[0] ) );
766  c10 = ( c010 * ( 1.0 - static_cast< proshade_double > ( interpDiff[0] ) ) ) + ( c110 * static_cast< proshade_double > ( interpDiff[0] ) );
767  c11 = ( c011 * ( 1.0 - static_cast< proshade_double > ( interpDiff[0] ) ) ) + ( c111 * static_cast< proshade_double > ( interpDiff[0] ) );
768 
769  //==================================== Interpolate along y-axis
770  c0 = ( c00 * ( 1.0 - static_cast< proshade_double > ( interpDiff[1] ) ) ) + ( c10 * static_cast< proshade_double > ( interpDiff[1] ) );
771  c1 = ( c01 * ( 1.0 - static_cast< proshade_double > ( interpDiff[1] ) ) ) + ( c11 * static_cast< proshade_double > ( interpDiff[1] ) );
772 
773  //==================================== Interpolate along z-axis
774  arrPos = static_cast< size_t > ( ( zIt - mins[2] ) + static_cast< proshade_single > ( this->zDimIndices ) * ( ( yIt - mins[1] ) + static_cast< proshade_single > ( this->yDimIndices ) * ( xIt - mins[0] ) ) );
775  map[arrPos] = ( c0 * ( 1.0 - static_cast< proshade_double > ( interpDiff[2] ) ) ) + ( c1 * static_cast< proshade_double > ( interpDiff[2] ) );
776  }
777  }
778  }
779 
780  //================================================ Release memory
781  delete[] mins;
782  delete[] maxs;
783  delete[] rotMat;
784  delete[] interpMins;
785  delete[] interpMaxs;
786  delete[] interpDiff;
787  delete[] movs;
788 
789  //================================================ Done
790  return ( ret );
791 
792 }
793 
805 std::vector< proshade_double > ProSHADE_internal_data::ProSHADE_data::rotateMapRealSpaceInPlace ( proshade_double eulA, proshade_double eulB, proshade_double eulG )
806 {
807  //================================================ Initialise local variables
808  proshade_double axX, axY, axZ, axAng, tmp, *rMat, *map;
809 
810  //================================================ Allocate local memory
811  rMat = new proshade_double[9];
812 
813  //================================================ Check local memory allocation
814  ProSHADE_internal_misc::checkMemoryAllocation ( rMat, __FILE__, __LINE__, __func__ );
815 
816  //================================================ Convert Euler angles to rotation matrix
818 
819  //================================================ Transpose the rotation matrix
820  tmp = rMat[1];
821  rMat[1] = rMat[3];
822  rMat[3] = tmp;
823 
824  tmp = rMat[2];
825  rMat[2] = rMat[6];
826  rMat[6] = tmp;
827 
828  tmp = rMat[5];
829  rMat[5] = rMat[7];
830  rMat[7] = tmp;
831 
832  //================================================ Convert rotation matrix to angle-axis
833  ProSHADE_internal_maths::getAxisAngleFromRotationMatrix ( rMat, &axX, &axY, &axZ, &axAng );
834 
835  //================================================ Rotate the internal map
836  std::vector< proshade_double > ret = this->rotateMapRealSpace ( axX, axY, axZ, axAng, map );
837 
838  //================================================ Copy the rotated map in place of the internal map
839  for ( size_t iter = 0; iter < static_cast< size_t > ( this->xDimIndices * this->yDimIndices * this->zDimIndices ); iter++ )
840  {
841  this->internalMap[iter] = map[iter];
842  }
843 
844  //================================================ Release local memory
845  delete[] rMat;
846  delete[] map;
847 
848  //================================================ Save rotation centre for co-ordinates
849  this->originalPdbRotCenX = ret.at(0);
850  this->originalPdbRotCenY = ret.at(1);
851  this->originalPdbRotCenZ = ret.at(2);
852 
853  //================================================ Done
854  return ( ret );
855 
856 }
857 
868 void ProSHADE_internal_data::ProSHADE_data::translateMap ( proshade_double trsX, proshade_double trsY, proshade_double trsZ )
869 {
870  //================================================ Initialise local variables
871  proshade_single xMov = static_cast< proshade_single > ( -trsX );
872  proshade_single yMov = static_cast< proshade_single > ( -trsY );
873  proshade_single zMov = static_cast< proshade_single > ( -trsZ );
874 
875  //================================================ Move the whole map frame to minimise the Fourier movement
876  ProSHADE_internal_mapManip::moveMapByIndices ( &xMov, &yMov, &zMov, this->getXDimSize(), this->getYDimSize(), this->getZDimSize(),
877  this->getXFromPtr(), this->getXToPtr(), this->getYFromPtr(), this->getYToPtr(),
878  this->getZFromPtr(), this->getZToPtr(), this->getXAxisOrigin(), this->getYAxisOrigin(), this->getZAxisOrigin() );
879 
880  //================================================ Finalise the movement by in-frame Fourier movement
881  ProSHADE_internal_mapManip::moveMapByFourier ( this->getInternalMap(), xMov, yMov, zMov, this->getXDimSize(), this->getYDimSize(), this->getZDimSize(),
882  static_cast< proshade_signed > ( this->getXDim() ), static_cast< proshade_signed > ( this->getYDim() ),
883  static_cast< proshade_signed > ( this->getZDim() ) );
884 
885  //================================================ Done
886  return ;
887 
888 }
889 
893 {
894  //================================================ Allocate the main pointer and check
895  this->rotSphericalHarmonics = new proshade_complex* [this->noSpheres];
896  ProSHADE_internal_misc::checkMemoryAllocation ( this->rotSphericalHarmonics, __FILE__, __LINE__, __func__ );
897 
898  //================================================ For each sphere
899  for ( proshade_unsign iter = 0; iter < this->noSpheres; iter++ )
900  {
901  //============================================ Allocate the sphere storage space
902  this->rotSphericalHarmonics[iter] = new proshade_complex [static_cast<proshade_unsign> ( pow ( (this->spheres[iter]->getLocalBandwidth() * 2), 2) )];
903  ProSHADE_internal_misc::checkMemoryAllocation ( this->rotSphericalHarmonics[iter], __FILE__, __LINE__, __func__ );
904 
905  //============================================ Set values to zeroes
906  for ( proshade_unsign it = 0; it < static_cast<proshade_unsign> ( pow ( (this->spheres[iter]->getLocalBandwidth() * 2), 2) ); it++ )
907  {
908  this->rotSphericalHarmonics[iter][it][0] = 0.0;
909  this->rotSphericalHarmonics[iter][it][1] = 0.0;
910  }
911  }
912 
913  //================================================ Done
914  return ;
915 
916 }
917 
921 {
922  //================================================ Initialise variables
923  proshade_double WigDR, WigDI, *ShR, *ShI, retR, retI;
924  proshade_unsign arrPos;
925 
926  //================================================ Compute
927  for ( proshade_signed shell = 0; shell < static_cast<proshade_signed> ( this->noSpheres ); shell++ )
928  {
929  //============================================ For each band
930  for ( proshade_signed bandIter = 0; bandIter < static_cast<proshade_signed> ( this->spheres[shell]->getLocalBandwidth() ); bandIter++ )
931  {
932  //======================================== For each order1
933  for ( proshade_signed order1 = 0; order1 < ( ( bandIter * 2 ) + 1 ); order1++ )
934  {
935  //==================================== Get Spherical Harmonics value
936  ShR = getRealSphHarmValue ( static_cast< proshade_unsign > ( bandIter ), static_cast< proshade_unsign > ( order1 ), static_cast< proshade_unsign > ( shell ) );
937  ShI = getImagSphHarmValue ( static_cast< proshade_unsign > ( bandIter ), static_cast< proshade_unsign > ( order1 ), static_cast< proshade_unsign > ( shell ) );
938 
939  //==================================== For each order2
940  for ( proshade_signed order2 = 0; order2 < ( ( bandIter * 2 ) + 1 ); order2++ )
941  {
942  //================================ Get Wigner D values
943  this->getWignerMatrixValue ( static_cast< proshade_unsign > ( bandIter ), static_cast< proshade_unsign > ( order1 ), static_cast< proshade_unsign > ( order2 ), &WigDR, &WigDI );
944 
945  //================================ Multiply SH and Wigner
946  ProSHADE_internal_maths::complexMultiplication ( ShR, ShI, &WigDR, &WigDI, &retR, &retI );
947 
948  //================================ Save
949  arrPos = static_cast<proshade_unsign> ( seanindex ( static_cast< int > ( order2-bandIter ), static_cast< int > ( bandIter ),
950  static_cast< int > ( this->spheres[shell]->getLocalBandwidth() ) ) );
951  this->rotSphericalHarmonics[shell][arrPos][0] += retR;
952  this->rotSphericalHarmonics[shell][arrPos][1] += retI;
953  }
954  }
955  }
956  }
957 
958  //================================================ Done
959  return ;
960 
961 }
962 
975 void ProSHADE_internal_overlay::initialiseInverseSHComputation ( proshade_unsign shBand, double*& sigR, double*& sigI, double*& rcoeffs, double*& icoeffs, double*& weights, double*& workspace, fftw_plan& idctPlan, fftw_plan& ifftPlan )
976 {
977  //================================================ Initialise internal variables
978  proshade_unsign oneDim = shBand * 2;
979 
980  //================================================ Allocate memory
981  sigR = new double [(oneDim * oneDim * 4)];
982  sigI = sigR + (oneDim * oneDim * 2);
983  rcoeffs = new double [(oneDim * oneDim * 2)];
984  icoeffs = rcoeffs + (oneDim * oneDim);
985  weights = new double [4 * shBand];
986  workspace = new double [( 20 * shBand * shBand ) + ( 24 * shBand )];
987 
988  //================================================ Check memory allocation
989  ProSHADE_internal_misc::checkMemoryAllocation ( sigR, __FILE__, __LINE__, __func__ );
990  ProSHADE_internal_misc::checkMemoryAllocation ( rcoeffs, __FILE__, __LINE__, __func__ );
991  ProSHADE_internal_misc::checkMemoryAllocation ( weights, __FILE__, __LINE__, __func__ );
992  ProSHADE_internal_misc::checkMemoryAllocation ( workspace, __FILE__, __LINE__, __func__ );
993 
994  //================================================ Create the cosine/sine transform plan
995  idctPlan = fftw_plan_r2r_1d ( static_cast< int > ( oneDim ), weights, workspace, FFTW_REDFT01, FFTW_ESTIMATE );
996 
997  //================================================ Set up the discrete Fourier transform
998  int rank, howmany_rank;
999  fftw_iodim dims[1], howmany_dims[1];
1000 
1001  rank = 1;
1002  dims[0].n = 2 * static_cast< int > ( shBand );
1003  dims[0].is = 2 * static_cast< int > ( shBand );
1004  dims[0].os = 1;
1005  howmany_rank = 1;
1006  howmany_dims[0].n = 2 * static_cast< int > ( shBand );
1007  howmany_dims[0].is = 1;
1008  howmany_dims[0].os = 2 * static_cast< int > ( shBand );
1009 
1010  //================================================ Create the discrete Fourier transform
1011  ifftPlan = fftw_plan_guru_split_dft ( rank, dims, howmany_rank, howmany_dims, sigR, sigI, rcoeffs, icoeffs, FFTW_ESTIMATE );
1012 
1013  //================================================ Done
1014  return ;
1015 
1016 }
1017 
1022 {
1023  //================================================ Initialise local variables
1024  double *sigR = nullptr, *sigI = nullptr, *rcoeffs = nullptr, *icoeffs = nullptr, *weights = nullptr, *workspace = nullptr;
1025  fftw_plan idctPlan, ifftPlan;
1026 
1027  //================================================ For each shell
1028  for ( int shell = 0; shell < static_cast<int> ( this->noSpheres ); shell++ )
1029  {
1030  //=========================================== Initialise internal variables
1031  proshade_unsign oneDim = this->spheres[shell]->getLocalBandwidth() * 2;
1032 
1033  //=========================================== Allocate memory
1034  ProSHADE_internal_overlay::initialiseInverseSHComputation ( this->spheres[shell]->getLocalBandwidth(), sigR, sigI, rcoeffs, icoeffs, weights, workspace, idctPlan, ifftPlan );
1035 
1036  //=========================================== Compute weights for the transform using the appropriate shell related band
1037  makeweights ( static_cast< int > ( this->spheres[shell]->getLocalBandwidth() ), weights );
1038 
1039  //============================================ Allocate rotated shell mapped data memory
1040  this->spheres[shell]->allocateRotatedMap ( );
1041 
1042  //============================================ Load SH coeffs to arrays
1043  for ( unsigned int iter = 0; iter < static_cast<unsigned int> ( oneDim * oneDim ); iter++ )
1044  {
1045  rcoeffs[iter] = this->rotSphericalHarmonics[shell][iter][0];
1046  icoeffs[iter] = this->rotSphericalHarmonics[shell][iter][1];
1047  sigR[iter] = 0.0;
1048  sigI[iter] = 0.0;
1049  }
1050 
1051  //============================================ Get inverse spherical harmonics transform for the shell
1052  InvFST_semi_fly ( rcoeffs,
1053  icoeffs,
1054  sigR,
1055  sigI,
1056  static_cast< int > ( this->spheres[shell]->getLocalBandwidth() ),
1057  workspace,
1058  0,
1059  static_cast< int > ( this->spheres[shell]->getLocalBandwidth() ),
1060  &idctPlan,
1061  &ifftPlan );
1062 
1063  //=========================================== Copy the results to the rotated shells array
1064  for ( unsigned int iter = 0; iter < static_cast<unsigned int> ( oneDim * oneDim ); iter++ )
1065  {
1066  this->spheres[shell]->setRotatedMappedData ( iter, sigR[iter] );
1067  }
1068 
1069  //=========================================== Release the plans
1070  fftw_destroy_plan ( idctPlan );
1071  fftw_destroy_plan ( ifftPlan );
1072 
1073  //=========================================== Release the memory
1074  delete[] sigR;
1075  delete[] rcoeffs;
1076  delete[] weights;
1077  delete[] workspace;
1078  }
1079 
1080  //================================================ Done
1081  return ;
1082 
1083 }
1084 
1091 void ProSHADE_internal_overlay::computeAngularThreshold ( std::vector<proshade_double>* lonCO, std::vector<proshade_double>* latCO, proshade_unsign angRes )
1092 {
1093  //================================================ Longitude angular thresholds
1094  for ( proshade_unsign iter = 0; iter <= angRes; iter++ )
1095  {
1096  ProSHADE_internal_misc::addToDoubleVector ( lonCO, static_cast<proshade_double> ( iter ) * ( ( static_cast<proshade_double> ( M_PI ) * 2.0 ) / static_cast<proshade_double> ( angRes ) ) - ( static_cast<double> ( M_PI ) ) );
1097  }
1098 
1099  //================================================ Lattitude angular thresholds
1100  for ( proshade_unsign iter = 0; iter <= angRes; iter++ )
1101  {
1102  ProSHADE_internal_misc::addToDoubleVector ( latCO, static_cast<proshade_double> ( iter ) * ( static_cast<proshade_double> ( M_PI ) / static_cast<proshade_double> ( angRes ) ) - ( static_cast<proshade_double> ( M_PI ) / 2.0 ) );
1103  }
1104 
1105  //================================================ Done
1106  return ;
1107 
1108 }
1109 
1114 void ProSHADE_internal_data::ProSHADE_data::interpolateMapFromSpheres ( proshade_double*& densityMapRotated )
1115 {
1116  //================================================ Initialise variables
1117  proshade_double rad = 0.0, lon = 0.0, lat = 0.0, newU = 0.0, newV = 0.0, newW = 0.0;
1118  proshade_unsign lowerLonL = 0, upperLonL = 0, lowerLonU = 0, upperLonU = 0, lowerLatL = 0, upperLatL = 0, lowerLatU = 0, upperLatU = 0, lowerShell = 0, upperShell = 0;
1119  proshade_double x00 = 0.0, x01 = 0.0, x10 = 0.0, x11 = 0.0, distLLon = 0.0, distLLat = 0.0, distLRad = 0.0, valLLon = 0.0, valULon = 0.0;
1120  proshade_double lowerShellValue = 0.0, upperShellValue = 0.0;
1121  proshade_double xSamplingRate = static_cast<proshade_double> ( this->xDimSize ) / static_cast<proshade_double> ( this->xDimIndices );
1122  proshade_double ySamplingRate = static_cast<proshade_double> ( this->yDimSize ) / static_cast<proshade_double> ( this->yDimIndices );
1123  proshade_double zSamplingRate = static_cast<proshade_double> ( this->zDimSize ) / static_cast<proshade_double> ( this->zDimIndices );
1124  proshade_signed arrPos;
1125  std::vector<proshade_double> lonCOU, latCOU, lonCOL, latCOL;
1126 
1127  for ( proshade_signed uIt = 0; uIt < static_cast<proshade_signed> (this->xDimIndices); uIt++ )
1128  {
1129  for ( proshade_signed vIt = 0; vIt < static_cast<proshade_signed> (this->yDimIndices); vIt++ )
1130  {
1131  for ( proshade_signed wIt = 0; wIt < static_cast<proshade_signed> (this->zDimIndices); wIt++ )
1132  {
1133  //==================================== Convert to centered coords
1134  newU = static_cast<proshade_double> ( uIt - ( static_cast<proshade_signed> (this->xDimIndices) / 2 ) );
1135  newV = static_cast<proshade_double> ( vIt - ( static_cast<proshade_signed> (this->yDimIndices) / 2 ) );
1136  newW = static_cast<proshade_double> ( wIt - ( static_cast<proshade_signed> (this->zDimIndices) / 2 ) );
1137 
1138  //==================================== Deal with 0 ; 0 ; 0
1139  if ( ( newU == 0.0 ) && ( newV == 0.0 ) && ( newW == 0.0 ) )
1140  {
1141  arrPos = wIt + static_cast< proshade_signed > ( this->zDimIndices ) * ( vIt + static_cast< proshade_signed > ( this->yDimIndices ) * uIt );
1142  densityMapRotated[arrPos] = this->internalMap[arrPos];
1143  continue;
1144  }
1145 
1146  //==================================== Convert to spherical coords
1147  rad = sqrt ( pow( ( newU * xSamplingRate ), 2.0 ) +
1148  pow( ( newV * ySamplingRate ), 2.0 ) +
1149  pow( ( newW * zSamplingRate ), 2.0 ) );
1150  lon = atan2 ( ( newV * ySamplingRate ), ( newU * xSamplingRate ) );
1151  lat = asin ( ( newW * zSamplingRate ) / rad );
1152 
1153  //==================================== Deal with nan's
1154  if ( rad != rad ) { rad = 0.0; }
1155  if ( lon != lon ) { lon = 0.0; }
1156  if ( lat != lat ) { lat = 0.0; }
1157 
1158  //==================================== Find shells above and below
1159  lowerShell = 0;
1160  upperShell = 0;
1161  for ( proshade_unsign iter = 0; iter < (this->noSpheres-1); iter++ )
1162  {
1163  if ( ( static_cast< proshade_double > ( this->spherePos.at(iter) ) <= rad ) && ( static_cast< proshade_double > ( this->spherePos.at(iter+1) ) > rad ) )
1164  {
1165  lowerShell = iter;
1166  upperShell = iter+1;
1167  break;
1168  }
1169  }
1170 
1171  if ( upperShell == 0 )
1172  {
1173  arrPos = wIt + static_cast< proshade_signed > ( this->zDimIndices ) * ( vIt + static_cast< proshade_signed > ( this->yDimIndices ) * uIt );
1174  densityMapRotated[arrPos] = 0.0;
1175  continue;
1176  }
1177 
1178  //==================================== Get the longitude and lattitude cut-offs for this shell resolution
1179  lonCOL.clear(); latCOL.clear(); lonCOU.clear(); latCOU.clear();
1180  ProSHADE_internal_overlay::computeAngularThreshold ( &lonCOL, &latCOL, this->spheres[lowerShell]->getLocalAngRes() );
1181  ProSHADE_internal_overlay::computeAngularThreshold ( &lonCOU, &latCOU, this->spheres[upperShell]->getLocalAngRes() );
1182 
1183  //==================================== Find the angle cutoffs around the point
1184  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( lonCOL.size() ); iter++ )
1185  {
1186  if ( iter == ( static_cast<proshade_unsign> ( lonCOL.size() ) - 1 ) )
1187  {
1188  lowerLonL = 0;
1189  upperLonL = 1;
1190  break;
1191  }
1192  if ( ( std::floor(10000. * lonCOL.at(iter)) <= std::floor(10000. * lon) ) && ( std::floor(10000. * lonCOL.at(iter+1)) > std::floor(10000. * lon) ) )
1193  {
1194  lowerLonL = iter;
1195  upperLonL = iter+1;
1196  break;
1197  }
1198  }
1199  if ( upperLonL == this->spheres[lowerShell]->getLocalAngRes() ) { upperLonL = 0; }
1200 
1201  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( lonCOU.size() ); iter++ )
1202  {
1203  if ( iter == ( static_cast<proshade_unsign> ( lonCOU.size() ) - 1 ) )
1204  {
1205  lowerLonU = 0;
1206  upperLonU = 1;
1207  break;
1208  }
1209  if ( ( std::floor(10000. * lonCOU.at(iter)) <= std::floor(10000. * lon) ) && ( std::floor(10000. * lonCOU.at(iter+1)) > std::floor(10000. * lon) ) )
1210  {
1211  lowerLonU = iter;
1212  upperLonU = iter+1;
1213  break;
1214  }
1215  }
1216  if ( upperLonU == this->spheres[upperShell]->getLocalAngRes() ) { upperLonU = 0; }
1217 
1218  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( latCOL.size() ); iter++ )
1219  {
1220  if ( iter == ( static_cast<proshade_unsign> ( latCOL.size() ) - 1 ) )
1221  {
1222  lowerLatL = 0;
1223  upperLatL = 1;
1224  break;
1225  }
1226  if ( ( std::floor(10000. * latCOL.at(iter)) <= std::floor(10000. * lat) ) && ( std::floor(10000. * latCOL.at(iter+1)) > std::floor(10000. * lat) ) )
1227  {
1228  lowerLatL = iter;
1229  upperLatL = iter+1;
1230  break;
1231  }
1232  }
1233  if ( upperLatL == this->spheres[lowerShell]->getLocalAngRes() ) { upperLatL = 0; }
1234 
1235  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( latCOU.size() ); iter++ )
1236  {
1237  if ( iter == ( static_cast<proshade_unsign> ( latCOU.size() ) - 1 ) )
1238  {
1239  lowerLatU = 0;
1240  upperLatU = 1;
1241  break;
1242  }
1243  if ( ( std::floor(10000. * latCOU.at(iter)) <= std::floor(10000. * lat) ) && ( std::floor(10000. * latCOU.at(iter+1)) > std::floor(10000. * lat) ) )
1244  {
1245  lowerLatU = iter;
1246  upperLatU = iter+1;
1247  break;
1248  }
1249  }
1250  if ( upperLatU == this->spheres[upperShell]->getLocalAngRes() ) { upperLatU = 0; }
1251 
1252  //==================================== Interpolate lower shell
1253  x00 = this->spheres[lowerShell]->getRotatedMappedData ( lowerLatL * this->spheres[lowerShell]->getLocalAngRes() + lowerLonL );
1254  x01 = this->spheres[lowerShell]->getRotatedMappedData ( lowerLatL * this->spheres[lowerShell]->getLocalAngRes() + upperLonL );
1255  x10 = this->spheres[lowerShell]->getRotatedMappedData ( upperLatL * this->spheres[lowerShell]->getLocalAngRes() + lowerLonL );
1256  x11 = this->spheres[lowerShell]->getRotatedMappedData ( upperLatL * this->spheres[lowerShell]->getLocalAngRes() + upperLonL );
1257 
1258  distLLon = std::abs ( lon - lonCOL.at(lowerLonL) ) / ( std::abs( lon - lonCOL.at(lowerLonL) ) + std::abs( lon - lonCOL.at(upperLonL) ) );
1259  valLLon = ( ( 1.0 - distLLon ) * x00 ) + ( distLLon * x01 );
1260  valULon = ( ( 1.0 - distLLon ) * x10 ) + ( distLLon * x11 );
1261 
1262  distLLat = std::abs ( lat - latCOL.at(lowerLatL) ) / ( std::abs( lat - latCOL.at(lowerLatL) ) + std::abs( lat - latCOL.at(upperLatL) ) );
1263  lowerShellValue = ( ( 1.0 - distLLat ) * valLLon ) + ( distLLat * valULon );
1264 
1265  //==================================== Interpolate upper shell
1266  x00 = this->spheres[upperShell]->getRotatedMappedData ( lowerLatU * this->spheres[upperShell]->getLocalAngRes() + lowerLonU );
1267  x01 = this->spheres[upperShell]->getRotatedMappedData ( lowerLatU * this->spheres[upperShell]->getLocalAngRes() + upperLonU );
1268  x10 = this->spheres[upperShell]->getRotatedMappedData ( upperLatU * this->spheres[upperShell]->getLocalAngRes() + lowerLonU );
1269  x11 = this->spheres[upperShell]->getRotatedMappedData ( upperLatU * this->spheres[upperShell]->getLocalAngRes() + upperLonU );
1270 
1271  distLLon = std::abs ( lon - lonCOU.at(lowerLonU) ) / ( std::abs( lon - lonCOU.at(lowerLonU) ) + std::abs( lon - lonCOU.at(upperLonU) ) );
1272  valLLon = ( ( 1.0 - distLLon ) * x00 ) + ( distLLon * x01 );
1273  valULon = ( ( 1.0 - distLLon ) * x10 ) + ( distLLon * x11 );
1274 
1275  distLLat = std::abs ( lat - latCOU.at(lowerLatU) ) / ( std::abs( lat - latCOU.at(lowerLatU) ) + std::abs( lat - latCOU.at(upperLatU) ) );
1276  upperShellValue = ( ( 1.0 - distLLat ) * valLLon ) + ( distLLat * valULon );
1277 
1278  //==================================== Interpolate between shells
1279  distLRad = std::abs ( rad - static_cast< proshade_double > ( this->spherePos.at(lowerShell) ) ) / ( std::abs( rad - static_cast< proshade_double > ( this->spherePos.at(lowerShell) ) ) +
1280  std::abs( rad - static_cast< proshade_double > ( this->spherePos.at(upperShell) ) ) );
1281 
1282  arrPos = wIt + static_cast< proshade_signed > ( this->zDimIndices ) * ( vIt + static_cast< proshade_signed > ( this->yDimIndices ) * uIt );
1283  densityMapRotated[arrPos] = ( ( 1.0 - distLRad ) * lowerShellValue ) + ( distLRad * upperShellValue );
1284  }
1285 
1286  }
1287 
1288  }
1289 
1290  //================================================ Done
1291  return ;
1292 
1293 }
1294 
1301 {
1302  //================================================ Initialise local variables
1303  std::vector < proshade_double > ret;
1304  proshade_double eulA, eulB, eulG;
1305 
1306  //================================================ Get inverse SO(3) map top peak Euler angle values
1308  this->getMaxBand() * 2,
1309  &eulA, &eulB, &eulG, settings );
1310 
1311  //================================================ Re-format to vector
1315 
1316  //================================================ Done
1317  return ( ret );
1318 
1319 }
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:572
ProSHADE_settings::maxBandwidth
proshade_unsign maxBandwidth
The bandwidth of spherical harmonics decomposition for the largest sphere.
Definition: ProSHADE_settings.hpp:58
ProSHADE_internal_data::ProSHADE_data::computeSphericalHarmonics
void computeSphericalHarmonics(ProSHADE_settings *settings)
This function computes the spherical harmonics decomposition for the whole structure.
Definition: ProSHADE_data.cpp:1847
ProSHADE_internal_data::ProSHADE_data::mapMovFromsChangeX
proshade_double mapMovFromsChangeX
When the map is translated, the xFrom and xTo values are changed. This variable holds how much they h...
Definition: ProSHADE_data.hpp:94
ProSHADE_internal_data::ProSHADE_data::allocateRotatedSHMemory
void allocateRotatedSHMemory(void)
This function allocates the memory required for storing the rotated Spherical Harmonics coefficients.
Definition: ProSHADE_overlay.cpp:892
ProSHADE_internal_data::ProSHADE_data::getOverlayRotationFunction
void getOverlayRotationFunction(ProSHADE_settings *settings, ProSHADE_internal_data::ProSHADE_data *obj2)
This function computes the overlay rotation function (i.e. the correlation function in SO(3) space).
Definition: ProSHADE_overlay.cpp:35
ProSHADE_internal_data::ProSHADE_data::zFrom
proshade_signed zFrom
This is the starting index along the z axis.
Definition: ProSHADE_data.hpp:112
ProSHADE_internal_data::ProSHADE_data::xDimIndices
proshade_unsign xDimIndices
This is the size of the map cell x dimension in indices.
Definition: ProSHADE_data.hpp:65
ProSHADE_internal_data::ProSHADE_data::getZDimSize
proshade_single getZDimSize(void)
This function allows access to the map size in angstroms along the Z axis.
Definition: ProSHADE_data.cpp:3919
ProSHADE_internal_data::ProSHADE_data::getBestTranslationMapPeaksAngstrom
std::vector< proshade_double > getBestTranslationMapPeaksAngstrom(ProSHADE_internal_data::ProSHADE_data *staticStructure)
This function gets the optimal translation vector and returns it as a standard library vector....
Definition: ProSHADE_overlay.cpp:250
ProSHADE_settings::determineAllSHValues
void determineAllSHValues(proshade_unsign xDim, proshade_unsign yDim, proshade_single xDimAngs, proshade_single yDimAngs, proshade_single zDimAngs)
This function determines all the required values for spherical harmonics computation.
Definition: ProSHADE.cpp:1853
ProSHADE_internal_data::ProSHADE_data::interpolateMapFromSpheres
void interpolateMapFromSpheres(proshade_double *&densityMapRotated)
This function interpolates the density map from the sphere mapped data.
Definition: ProSHADE_overlay.cpp:1114
ProSHADE_internal_data::ProSHADE_data::getXAxisOrigin
proshade_signed * getXAxisOrigin(void)
This function allows access to the map X axis origin value.
Definition: ProSHADE_data.cpp:4019
ProSHADE_exception
This class is the representation of ProSHADE exception.
Definition: ProSHADE_exceptions.hpp:37
ProSHADE_internal_data::ProSHADE_data::getXDimSize
proshade_single getXDimSize(void)
This function allows access to the map size in angstroms along the X axis.
Definition: ProSHADE_data.cpp:3899
ProSHADE_internal_data::ProSHADE_data::getInternalMap
proshade_double *& getInternalMap(void)
This function allows access to the first map array value address.
Definition: ProSHADE_data.cpp:4049
ProSHADE_internal_data::ProSHADE_data
This class contains all inputed and derived data for a single structure.
Definition: ProSHADE_data.hpp:49
ProSHADE_internal_data::ProSHADE_data::getYDimSize
proshade_single getYDimSize(void)
This function allows access to the map size in angstroms along the Y axis.
Definition: ProSHADE_data.cpp:3909
ProSHADE_internal_data::ProSHADE_data::zDimSize
proshade_single zDimSize
This is the size of the map cell z dimension in Angstroms.
Definition: ProSHADE_data.hpp:61
ProSHADE_internal_mapManip::findMAPCOMValues
void findMAPCOMValues(proshade_double *map, proshade_double *xCom, proshade_double *yCom, proshade_double *zCom, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_signed xFrom, proshade_signed xTo, proshade_signed yFrom, proshade_signed yTo, proshade_signed zFrom, proshade_signed zTo)
This function finds the Centre of Mass for a map.
Definition: ProSHADE_mapManip.cpp:240
ProSHADE_internal_data::ProSHADE_data::getBestRotationMapPeaksEulerAngles
std::vector< proshade_double > getBestRotationMapPeaksEulerAngles(ProSHADE_settings *settings)
This function returns a vector of three floats, the three Euler angles of the best peak in the rotati...
Definition: ProSHADE_overlay.cpp:1300
ProSHADE_internal_mapManip::moveMapByFourier
void moveMapByFourier(proshade_double *&map, proshade_single xMov, proshade_single yMov, proshade_single zMov, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_signed xDim, proshade_signed yDim, proshade_signed zDim)
Function for moving map back to original PDB location by using Fourier transformation.
Definition: ProSHADE_mapManip.cpp:813
ProSHADE_internal_maths::getRotationMatrixFromEulerZXZAngles
void getRotationMatrixFromEulerZXZAngles(proshade_double eulerAlpha, proshade_double eulerBeta, proshade_double eulerGamma, proshade_double *matrix)
Function to find the rotation matrix from Euler angles (ZXZ convention).
Definition: ProSHADE_maths.cpp:1020
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:492
ProSHADE_internal_maths::complexMultiplication
void complexMultiplication(proshade_double *r1, proshade_double *i1, proshade_double *r2, proshade_double *i2, proshade_double *retReal, proshade_double *retImag)
Function to multiply two complex numbers.
Definition: ProSHADE_maths.cpp:38
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:74
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:406
ProSHADE_internal_data::ProSHADE_data::getXFromPtr
proshade_signed * getXFromPtr(void)
This function allows access to the map start along the X axis.
Definition: ProSHADE_data.cpp:3959
ProSHADE_settings::verbose
proshade_signed verbose
Should the software report on the progress, or just be quiet? Value between -1 (nothing) and 4 (loud)
Definition: ProSHADE_settings.hpp:149
ProSHADE_internal_data::ProSHADE_data::yDimIndices
proshade_unsign yDimIndices
This is the size of the map cell y dimension in indices.
Definition: ProSHADE_data.hpp:66
ProSHADE_settings::messageShift
proshade_signed messageShift
This value allows shifting the messages to create more readable log for sub-processes.
Definition: ProSHADE_settings.hpp:150
ProSHADE_internal_misc::addToDoubleVector
void addToDoubleVector(std::vector< proshade_double > *vecToAddTo, proshade_double elementToAdd)
Adds the element to the vector.
Definition: ProSHADE_misc.cpp:77
ProSHADE_internal_data::ProSHADE_data::getXDim
proshade_unsign getXDim(void)
This function allows access to the map size in indices along the X axis.
Definition: ProSHADE_data.cpp:3929
ProSHADE_internal_data::ProSHADE_data::originalPdbTransX
proshade_double originalPdbTransX
The optimal translation vector as it relates to the original PDB positions (and not the ProSHADE inte...
Definition: ProSHADE_data.hpp:105
ProSHADE_internal_data::ProSHADE_data::rotateMapReciprocalSpace
void rotateMapReciprocalSpace(ProSHADE_settings *settings, proshade_double eulerAlpha, proshade_double eulerBeta, proshade_double eulerGamma)
This function rotates a map based on the given Euler angles.
Definition: ProSHADE_overlay.cpp:571
ProSHADE_internal_data::ProSHADE_data::readInStructure
void readInStructure(std::string fName, proshade_unsign inputO, ProSHADE_settings *settings, proshade_double *maskArr=nullptr, proshade_unsign maskXDim=0, proshade_unsign maskYDim=0, proshade_unsign maskZDim=0, proshade_double *weightsArr=nullptr, proshade_unsign weigXDim=0, proshade_unsign weigYDim=0, proshade_unsign weigZDim=0)
This function initialises the basic ProSHADE_data variables and reads in a single structure.
Definition: ProSHADE_data.cpp:509
ProSHADE_internal_data::ProSHADE_data::mapToSpheres
void mapToSpheres(ProSHADE_settings *settings)
This function converts the internal map onto a set of concentric spheres.
Definition: ProSHADE_data.cpp:1803
ProSHADE_internal_data::ProSHADE_data::originalPdbTransY
proshade_double originalPdbTransY
The optimal translation vector as it relates to the original PDB positions (and not the ProSHADE inte...
Definition: ProSHADE_data.hpp:106
ProSHADE_internal_wigner::computeWignerMatricesForRotation
void computeWignerMatricesForRotation(ProSHADE_settings *settings, ProSHADE_internal_data::ProSHADE_data *obj, proshade_double eulerAlpha, proshade_double eulerBeta, proshade_double eulerGamma)
This function computes the Wigner D matrices for a particular set of Euler angles.
Definition: ProSHADE_wignerMatrices.cpp:258
ProSHADE_overlay.hpp
This header file declares the functions required for the structure overlay computation.
ProSHADE_internal_data::ProSHADE_data::xDimSize
proshade_single xDimSize
This is the size of the map cell x dimension in Angstroms.
Definition: ProSHADE_data.hpp:59
ProSHADE_internal_data::ProSHADE_data::computeRotatedSH
void computeRotatedSH(void)
This function multiplies the objects spherical harmonics with the Wigner D matrices,...
Definition: ProSHADE_overlay.cpp:920
ProSHADE_internal_data::ProSHADE_data::mapMovFromsChangeY
proshade_double mapMovFromsChangeY
When the map is translated, the yFrom and yTo values are changed. This variable holds how much they h...
Definition: ProSHADE_data.hpp:95
ProSHADE_internal_peakSearch::getBestPeakEulerAngsNaive
void getBestPeakEulerAngsNaive(proshade_complex *map, proshade_unsign dim, proshade_double *eulA, proshade_double *eulB, proshade_double *eulG, ProSHADE_settings *settings)
This function finds the highest peaks optimised Euler angles using the "naive" approach.
Definition: ProSHADE_peakSearch.cpp:352
ProSHADE_internal_data::ProSHADE_data::computeTranslationMap
void computeTranslationMap(ProSHADE_internal_data::ProSHADE_data *obj1)
This function does the computation of the translation map and saves results internally.
Definition: ProSHADE_overlay.cpp:324
ProSHADE_internal_data::ProSHADE_data::xFrom
proshade_signed xFrom
This is the starting index along the x axis.
Definition: ProSHADE_data.hpp:110
ProSHADE_internal_maths::getRotationMatrixFromAngleAxis
void getRotationMatrixFromAngleAxis(proshade_double *rotMat, proshade_double x, proshade_double y, proshade_double z, proshade_double ang)
This function converts the axis-angle representation to the rotation matrix representation.
Definition: ProSHADE_maths.cpp:1459
ProSHADE_internal_data::ProSHADE_data::getMaxBand
proshade_unsign getMaxBand(void)
This function returns the maximum band value for the object.
Definition: ProSHADE_data.cpp:3611
ProSHADE_internal_data::ProSHADE_data::translateMap
void translateMap(proshade_double trsX, proshade_double trsY, proshade_double trsZ)
This function simply translates the map by a given number of Angstroms along the three axes.
Definition: ProSHADE_overlay.cpp:868
ProSHADE_internal_data::ProSHADE_data::invertSHCoefficients
void invertSHCoefficients(void)
This function computes the shell mapped data from inverting the Spherical Harmonics coefficients.
Definition: ProSHADE_overlay.cpp:1021
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:123
ProSHADE_internal_data::ProSHADE_data::getYDim
proshade_unsign getYDim(void)
This function allows access to the map size in indices along the Y axis.
Definition: ProSHADE_data.cpp:3939
ProSHADE_internal_maths::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_maths.cpp:3912
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:37
ProSHADE_internal_maths::getAxisAngleFromRotationMatrix
void getAxisAngleFromRotationMatrix(proshade_double *rotMat, proshade_double *x, proshade_double *y, proshade_double *z, proshade_double *ang, proshade_signed verbose=1)
This function converts rotation matrix to the axis-angle representation.
Definition: ProSHADE_maths.cpp:1084
ProSHADE_internal_data::ProSHADE_data::originalPdbTransZ
proshade_double originalPdbTransZ
The optimal translation vector as it relates to the original PDB positions (and not the ProSHADE inte...
Definition: ProSHADE_data.hpp:107
ProSHADE_internal_data::ProSHADE_data::processInternalMap
void processInternalMap(ProSHADE_settings *settings)
This function simply clusters several other functions which should be called together.
Definition: ProSHADE_data.cpp:1695
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:510
ProSHADE_internal_data::ProSHADE_data::getTranslationFnPointer
proshade_complex * getTranslationFnPointer(void)
This function allows access to the translation function through a pointer.
Definition: ProSHADE_data.cpp:4059
ProSHADE_internal_maths::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_maths.cpp:3865
ProSHADE_internal_data::ProSHADE_data::rotateMapRealSpaceInPlace
std::vector< proshade_double > rotateMapRealSpaceInPlace(proshade_double eulA, proshade_double eulB, proshade_double eulG)
This function rotates a map based on the given Euler angles in place.
Definition: ProSHADE_overlay.cpp:805
ProSHADE_internal_data::ProSHADE_data::getZDim
proshade_unsign getZDim(void)
This function allows access to the map size in indices along the Z axis.
Definition: ProSHADE_data.cpp:3949
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:1091
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:975
ProSHADE_internal_data::ProSHADE_data::mapMovFromsChangeZ
proshade_double mapMovFromsChangeZ
When the map is translated, the zFrom and zTo values are changed. This variable holds how much they h...
Definition: ProSHADE_data.hpp:96
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:854
ProSHADE_internal_data::ProSHADE_data::getMapValue
proshade_double getMapValue(proshade_unsign pos)
This function returns the internal map representation value of a particular array position.
Definition: ProSHADE_data.cpp:3601
ProSHADE_internal_data::ProSHADE_data::getYToPtr
proshade_signed * getYToPtr(void)
This function allows access to the map last position along the Y axis.
Definition: ProSHADE_data.cpp:3999
ProSHADE_internal_data::ProSHADE_data::rotateMapRealSpace
std::vector< proshade_double > rotateMapRealSpace(proshade_double axX, proshade_double axY, proshade_double axZ, proshade_double axAng, proshade_double *&map)
This function rotates a map based on the given angle-axis rotation.
Definition: ProSHADE_overlay.cpp:632
ProSHADE_internal_data::ProSHADE_data::getZAxisOrigin
proshade_signed * getZAxisOrigin(void)
This function allows access to the map Z axis origin value.
Definition: ProSHADE_data.cpp:4039
ProSHADE_internal_maths::compute3x3MatrixVectorMultiplication
proshade_double * compute3x3MatrixVectorMultiplication(proshade_double *mat, proshade_double x, proshade_double y, proshade_double z)
Function for computing a 3x3 matrix to 3x1 vector multiplication.
Definition: ProSHADE_maths.cpp:1898
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:68
ProSHADE_internal_data::ProSHADE_data::getXToPtr
proshade_signed * getXToPtr(void)
This function allows access to the map last position along the X axis.
Definition: ProSHADE_data.cpp:3989
ProSHADE_internal_data::ProSHADE_data::getYFromPtr
proshade_signed * getYFromPtr(void)
This function allows access to the map start along the Y axis.
Definition: ProSHADE_data.cpp:3969
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:367
ProSHADE_internal_data::ProSHADE_data::getInvSO3Coeffs
proshade_complex * getInvSO3Coeffs(void)
This function allows access to the inverse SO(3) coefficients array.
Definition: ProSHADE_data.cpp:3847
ProSHADE_internal_data::ProSHADE_data::yFrom
proshade_signed yFrom
This is the starting index along the y axis.
Definition: ProSHADE_data.hpp:111
ProSHADE_settings::inputFiles
std::vector< std::string > inputFiles
This vector contains the filenames of all input structure files.
Definition: ProSHADE_settings.hpp:43
ProSHADE_internal_data::ProSHADE_data::getYAxisOrigin
proshade_signed * getYAxisOrigin(void)
This function allows access to the map Y axis origin value.
Definition: ProSHADE_data.cpp:4029
ProSHADE_internal_overlay::computeTranslationsFromPeak
void computeTranslationsFromPeak(ProSHADE_internal_data::ProSHADE_data *staticStructure, ProSHADE_internal_data::ProSHADE_data *movingStructure, proshade_double *trsX, proshade_double *trsY, proshade_double *trsZ)
This function computes the translation in Angstroms that corresponds to the translation function peak...
Definition: ProSHADE_overlay.cpp:220
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:521
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:701
ProSHADE_internal_messages::printProgressMessage
void printProgressMessage(proshade_signed verbose, proshade_signed messageLevel, std::string message, proshade_signed messageShift=0)
General stdout message printing.
Definition: ProSHADE_messages.cpp:71
ProSHADE_internal_data::ProSHADE_data::yDimSize
proshade_single yDimSize
This is the size of the map cell y dimension in Angstroms.
Definition: ProSHADE_data.hpp:60
ProSHADE_internal_mapManip::moveMapByIndices
void moveMapByIndices(proshade_single *xMov, proshade_single *yMov, proshade_single *zMov, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_signed *xFrom, proshade_signed *xTo, proshade_signed *yFrom, proshade_signed *yTo, proshade_signed *zFrom, proshade_signed *zTo, proshade_signed *xOrigin, proshade_signed *yOrigin, proshade_signed *zOrigin)
Function for moving map back to original PDB location by changing the indices.
Definition: ProSHADE_mapManip.cpp:765
ProSHADE_internal_data::ProSHADE_data::zeroPaddToDims
void zeroPaddToDims(proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim)
This function changes the size of a structure to fit the supplied new limits.
Definition: ProSHADE_overlay.cpp:433
ProSHADE_internal_data::ProSHADE_data::getZFromPtr
proshade_signed * getZFromPtr(void)
This function allows access to the map start along the Z axis.
Definition: ProSHADE_data.cpp:3979
ProSHADE_internal_data::ProSHADE_data::getZToPtr
proshade_signed * getZToPtr(void)
This function allows access to the map last position along the Z axis.
Definition: ProSHADE_data.cpp:4009