ProSHADE  0.6.5 (NOV 2018)
Protein Shape Descriptors and Symmetry Detection
ProSHADE_internal.h
Go to the documentation of this file.
1 
19 //============================================ ProSHADE
20 #include "ProSHADE_version.h"
21 #include "ProSHADE_rvapi.h"
22 
23 //============================================ Include once
24 #ifndef __PROSHADE_INTERNAL_LIBRARY__
25 #define __PROSHADE_INTERNAL_LIBRARY__
26 
33 namespace ProSHADE_internal
34 {
35  //======================================== Forward declaration
36  unsigned int checkFileType ( std::string fileName );
37 
38  //======================================== Forward declaration
39  class ProSHADE_comparePairwise;
40 
49  {
50  //==================================== Friends
51  friend class ProSHADE_comparePairwise;
52 
53  private:
54  //==================================== Settings regarding resolutions
55  double mapResolution;
56  std::vector<unsigned int> bandwidth;
57  std::vector<unsigned int> glIntegOrder;
58  std::vector<unsigned int> theta;
59  std::vector<unsigned int> phi;
60 
61  //==================================== Settings regarding B factors
62  double bFactorValue;
63  double bFactorChange;
64 
65  //==================================== Setting regarding maps and removing noise
66  double noIQRsFromMap;
67 
68  //==================================== Settings regarding concentric shells
69  double shellSpacing;
70  unsigned int manualShells;
71 
72  //==================================== Settings regarding map with phases
73  bool useCOM;
74  bool firstLineCOM;
75 
76  //==================================== Settings regarding space around the structure in lattice
77  double extraSpace;
78 
79  //==================================== Settings regarding weighting the distances
80  double alpha;
81  double mPower;
82 
83  //==================================== Settings regarding bands to ignore
84  std::vector<int> ignoreLs;
85 
86  //==================================== Settings regarding the structures to use
87  std::vector <std::string> structFiles;
88 
89  //==================================== Variables holding the results
90  std::vector< std::array<double,8> > rfPeaks;
91  std::vector< std::vector<std::array<double,5> > > tetrSymm;
92  double tetrSymmPeakAvg;
93  std::vector< std::array<double,5> > tetrElems;
94  std::vector< std::vector<std::array<double,5> > > octaSymm;
95  double octaSymmPeakAvg;
96  std::vector< std::array<double,5> > octaElems;
97  std::vector< std::vector<std::array<double,5> > > icosSymm;
98  double icosSymmPeakAvg;
99  std::vector< std::array<double,5> > icosElems;
100 
101  //==================================== Settings regarding peak finding
102  double peakHeightNoIQRs;
103  double peakDistanceForReal;
104  int peakSurroundingPoints;
105 
106  //==================================== Settings regarding axis and angle error tolerance
107  double aaErrorTolerance;
108  double symGapTolerance;
109 
110  //==================================== Are full results needed?
111  bool printFull;
112 
113  public:
114  //==================================== Public variables holding the results
115  std::vector< std::array<double,5> > cnSymm;
116  std::vector< std::array<double,5> > cnSymmClear;
117  std::vector< std::vector<std::array<double,6> > > dnSymm;
118  std::vector< std::vector<std::array<double,6> > > dnSymmClear;
119  std::vector< std::array<double,5> > tetrAxes;
120  std::vector< std::array<double,5> > octaAxes;
121  std::vector< std::array<double,5> > icosAxes;
123 
124  //==================================== Public functions
126  ProSHADE_symmetry ( );
127  void printResultsClear ( int verbose );
129  void printResultsRequest ( std::string symmetryType,
130  unsigned int symmetryFold,
131  int verbose );
132  void printResultsRequestHTML ( std::string symmetryType,
133  unsigned int symmetryFold,
134  ProSHADE::ProSHADE_settings* settings );
135  std::vector< std::array<double,8> > getRotFnPeaks ( );
136  std::vector< std::array<double,5> > getCSymmetries ( );
137  std::vector< std::vector<std::array<double,6> > > getDSymmetries ( );
138  std::vector< std::array<double,5> > generateTetrAxes ( ProSHADE_comparePairwise* cmpObj,
139  ProSHADE::ProSHADE_settings* settings,
140  std::vector< std::array<double,5> > tetrSymm,
141  std::vector< std::array<double,5> > allCs,
142  double axisErrorTolerance = 0.1,
143  int verbose = 0 );
144  std::vector< std::array<double,5> > generateTetrElements ( std::vector< std::array<double,5> > symmAxes,
145  ProSHADE::ProSHADE_settings* settings,
146  int verbose = 0 );
147  std::vector< std::array<double,5> > generateOctaAxes ( ProSHADE_comparePairwise* cmpObj,
148  ProSHADE::ProSHADE_settings* settings,
149  std::vector< std::array<double,5> > octaSymm,
150  std::vector< std::array<double,5> > allCs,
151  double axisErrorTolerance = 0.1,
152  int verbose = 0 );
153  std::vector< std::array<double,5> > generateOctaElements ( std::vector< std::array<double,5> > symmAxes,
154  ProSHADE::ProSHADE_settings* settings,
155  int verbose = 0 );
156  std::vector< std::array<double,5> > generateIcosAxes ( ProSHADE_comparePairwise* cmpObj,
157  ProSHADE::ProSHADE_settings* settings,
158  std::vector< std::array<double,5> > icosSymm,
159  std::vector< std::array<double,5> > allCs,
160  double axisErrorTolerance = 0.1,
161  int verbose = 0 );
162  std::vector< std::array<double,5> > generateIcosElements ( std::vector< std::array<double,5> > symmAxes,
163  ProSHADE::ProSHADE_settings* settings,
164  int verbose = 0 );
165  };
166 
175  {
176  //==================================== Friends
177  friend class ProSHADE_comparePairwise;
178  friend class ProSHADE_compareOneAgainstAll;
179  friend class ProSHADE_symmetry;
180  friend class ProSHADE_distances;
181  friend class ProSHADE_mapFeatures;
182 
183  private:
184  //==================================== Variables regarding the reading in of the data
185  std::string _inputFileName;
186  bool _fromPDB;
187 
188  //==================================== Variables regarding the shells
189  double _shellSpacing;
190  double _maxExtraCellularSpace;
191  double _xRange;
192  double _yRange;
193  double _zRange;
194  double _xSamplingRate;
195  double _ySamplingRate;
196  double _zSamplingRate;
197  float _xFrom;
198  float _yFrom;
199  float _zFrom;
200  float _xTo;
201  float _yTo;
202  float _zTo;
203  std::vector<double> _shellPlacement;
204 
205  //==================================== Variables regarding the density map
206  double _mapResolution;
207  unsigned int _maxMapU;
208  unsigned int _maxMapV;
209  unsigned int _maxMapW;
210  float *_densityMapMap;
211  bool _densityMapComputed;
212  double _mapMean;
213  double _mapSdev;
214 
215  //==================================== Variables regarding the phase removal from map
216  double _fourierCoeffPower;
217  double _bFactorChange;
218  double *_densityMapCor;
219  double _maxMapRange;
220  std::array<double,3> *_densityMapCorCoords;
221  bool _phaseRemoved;
222  bool _usePhase;
223  bool _keepOrRemove;
224 
225  //==================================== Variables regarding the sphere mapping of phaseless data
226  double _thetaAngle;
227  double _phiAngle;
228  unsigned int _noShellsWithData;
229  double **_shellMappedData;
230  int _xCorrection;
231  int _yCorrection;
232  int _zCorrection;
233  double _xCorrErr;
234  double _yCorrErr;
235  double _zCorrErr;
236  bool _sphereMapped;
237  bool _firstLineCOM;
238 
239  //==================================== Variables regarding the spherical harmonics decomposition
240  unsigned int _bandwidthLimit;
241  unsigned int _oneDimmension;
242  double **_realSHCoeffs;
243  double **_imagSHCoeffs;
244  double *_sphericalHarmonicsWeights;
245  double **_semiNaiveTable;
246  double *_semiNaiveTableSpace;
247  fftw_complex *_shWorkspace;
248  bool _sphericalCoefficientsComputed;
249  bool _wasBandwithGiven;
250  bool _wasThetaGiven;
251  bool _wasPhiGiven;
252  bool _wasGlInterGiven;
253 
254  //==================================== Variables regarding spherical harmonics inverse
255  double **_invRealData;
256  double **_invImagData;
257 
258  //==================================== Variables regarding the energy levels descriptor pre-calculation
259  double ***_rrpMatrices;
260  bool _rrpMatricesPrecomputed;
261 
262  public:
263  //==================================== Public functions
264  ProSHADE_data ( );
265  ProSHADE_data ( ProSHADE_data *copyFrom );
266  ProSHADE_data ( std::ifstream* dbFile, std::string fName, double volThreMin, double volThreMax, int verbose, ProSHADE::ProSHADE_settings* settings, bool saveWithAndWithout = false );
267  ~ProSHADE_data ( );
268 
269  void getDensityMapFromPDB ( std::string fileName, double *shellDistance, double resolution,
270  unsigned int *bandwidth, unsigned int *theta,
271  unsigned int *phi, unsigned int *glIntegOrder,
272  double *extraSpace, bool mapResDefault, ProSHADE::ProSHADE_settings* settings,
273  double Bfactor = 80.0, bool hpFirstLineCom = false, bool overlayDefaults = false );
274  void getDensityMapFromMAP ( std::string fileName, double *shellDistance, double resolution,
275  unsigned int *bandwidth, unsigned int *theta,
276  unsigned int *phi, unsigned int *glIntegOrder,
277  double *extraSpace, bool mapResDefault, bool rotDefaults, ProSHADE::ProSHADE_settings* settings,
278  bool overlayDefaults = false );
279  void maskMap ( int hlpU, int hlpV, int hlpW, double blurBy = 250.0, double maxMapIQR = 3.0, double extraMaskSpace = 3.0 );
280  void normaliseMap ( ProSHADE::ProSHADE_settings* settings );
281  std::array<double,4> getCOMandDist ( ProSHADE::ProSHADE_settings* settings );
282  bool removeIslands ( int hlpU, int hlpV, int hlpW, int verbose = 0, bool runAll = false );
283  void translateMap ( double xShift, double yShift, double zShift );
284  void shiftMap ( double xShift, double yShift, double zShift );
285  void matchMap ( ProSHADE_data *matchTo );
286  std::array<double,6> getDensityMapFromMAPRebox ( std::string fileName, double maxMapIQR, double extraCS, int verbose,
287  bool useCubicMaps, double maskBlurFactor, bool maskBlurFactorGiven,
288  ProSHADE::ProSHADE_settings* settings );
289  void getDensityMapFromMAPFeatures ( std::string fileName, double *minDensPreNorm, double *maxDensPreNorm,
290  double *minDensPostNorm, double *maxDensPostNorm, double *postNormMean,
291  double *postNormSdev, double *maskVolume, double* totalVolume, double *maskMean,
292  double *maskSdev, double *maskMin, double *maskMax, double *maskDensityRMS,
293  double *allDensityRMS, std::array<double,3> *origRanges, std::array<double,3> *origDims,
294  double maxMapIQR, double extraCS, int verbose, bool useCubicMaps, double maskBlurFactor, bool maskBlurFactorGiven,
295  bool reboxAtAll, ProSHADE::ProSHADE_settings* settings );
296  void removePhaseFromMap ( double alpha, double bFac, ProSHADE::ProSHADE_settings* settings );
297  void removePhaseFromMapOverlay ( double alpha, double bFac, ProSHADE::ProSHADE_settings* settings );
298  void keepPhaseInMap ( double alpha, double bFac,
299  unsigned int *bandwidth, unsigned int *theta,
300  unsigned int *phi, unsigned int *glIntegOrder, ProSHADE::ProSHADE_settings* settings,
301  bool useCom = true,
302  double maxMapIQR = 10.0, int verbose = 0,
303  bool clearMapData = true, bool rotDefaults = false, bool overlapDefaults = false,
304  double blurFactor = 500.0, bool maskBlurFactorGiven = false );
305  void mapPhaselessToSphere ( ProSHADE::ProSHADE_settings* settings, double theta, double phi, double shellSz, unsigned int manualShells = 0, bool keepInMemory = false, bool rotDefaults = false );
306  void getSphericalHarmonicsCoeffs ( unsigned int bandwidth, ProSHADE::ProSHADE_settings* settings );
307  void precomputeRotInvDescriptor ( ProSHADE::ProSHADE_settings* settings );
308  std::vector<ProSHADE_data*> fragmentMap ( ProSHADE::ProSHADE_settings* settings, bool userCOM );
309  std::vector< std::vector<int> > findMapIslands ( int hlpU, int hlpV, int hlpW, double *map, double threshold = 0.0 );
310  void writeMap ( std::string fileName,
311  double *map,
312  std::string axOrder = "xyz",
313  float xFrom = std::numeric_limits<float>::infinity(),
314  float yFrom = std::numeric_limits<float>::infinity(),
315  float zFrom = std::numeric_limits<float>::infinity(),
316  float xTo = std::numeric_limits<float>::infinity(),
317  float yTo = std::numeric_limits<float>::infinity(),
318  float zTo = std::numeric_limits<float>::infinity(),
319  float xRange = std::numeric_limits<float>::infinity(),
320  float yRange = std::numeric_limits<float>::infinity(),
321  float zRange = std::numeric_limits<float>::infinity() );
322  void writeMap ( std::string fileName,
323  float *map,
324  std::string axOrder = "xyz",
325  float xFrom = std::numeric_limits<float>::infinity(),
326  float yFrom = std::numeric_limits<float>::infinity(),
327  float zFrom = std::numeric_limits<float>::infinity(),
328  float xTo = std::numeric_limits<float>::infinity(),
329  float yTo = std::numeric_limits<float>::infinity(),
330  float zTo = std::numeric_limits<float>::infinity(),
331  float xRange = std::numeric_limits<float>::infinity(),
332  float yRange = std::numeric_limits<float>::infinity(),
333  float zRange = std::numeric_limits<float>::infinity() );
334  void writePDB ( std::string templateName, std::string outputName,
335  double rotEulA = 0.0,
336  double rotEulB = 0.0,
337  double rotEulG = 0.0,
338  double trsX = 0.0,
339  double trsY = 0.0,
340  double trsZ = 0.0 );
341  void deleteModel ( std::string modelPath );
342 
343  //==================================== Accessor functions
344  inline double getMapXRange ( void ) { return ( this->_xRange ); };
345  inline double getMapYRange ( void ) { return ( this->_yRange ); };
346  inline double getMapZRange ( void ) { return ( this->_zRange ); };
347  inline double* getMap ( void ) { return ( this->_densityMapCor ); };
348 
349  };
350 
359  {
360  private:
361  //==================================== Variables regarding sanity check
362  unsigned int _bandwidthLimit;
363  double _shellSpacing;
364  double _thetaAngle;
365  double _phiAngle;
366  double _matrixPowerWeight;
367  unsigned int _minShellsToUse;
368  std::vector<int> _lsToIgnore;
369  bool _bothRRPsPreComputed;
370  bool _keepOrRemove;
371 
372  //==================================== Variables regarding RotInv distance
373  double ***_obj1RRPs;
374  double ***_obj2RRPs;
375  double _distanceRotInv;
376  bool _rotInvComputed;
377 
378  //==================================== Variables regarding TrSigma distance pre-calculation
379  unsigned int _noShellsObj1;
380  unsigned int _noShellsObj2;
381  unsigned int _maxShellsToUse;
382  double **_obj1RealCoeffs;
383  double **_obj1ImagCoeffs;
384  double **_obj2RealCoeffs;
385  double **_obj2ImagCoeffs;
386  std::vector<double> _trSigmaWeights;
387  std::array<double,2> ***_trSigmaEMatrix;
388  bool _trSigmaPreComputed;
389 
390  //==================================== Variables regarding TrSigma distance computation
391  double _distanceTrSigma;
392  bool _trSigmaComputed;
393 
394  //==================================== Variables regarding the SO3 Transform
395  fftw_complex *_so3Coeffs;
396  fftw_complex *_so3InvCoeffs;
397  fftw_complex *_so3Workspace1;
398  fftw_complex *_so3Workspace2;
399  double *_so3Workspace3;
400  bool _so3InvMapComputed;
401 
402  //==================================== Variables regarding the Euler angles
403  std::array<double,3> _eulerAngles;
404  bool _eulerAnglesFound;
405 
406  //==================================== Variables regarding symmetry detection
407  double _peakHeightThr;
408  bool _CSymmsFound;
409  bool _DSymmsFound;
410 
411  //==================================== Variables regarding the Wigner d matrices computation
412  std::vector< std::vector< std::vector< std::array<double,2> > > > _wignerMatrices;
413  bool _wignerMatricesComputed;
414 
415  //==================================== Variables regarding Gauss-Legendre integration
416  unsigned int _glIntegrationOrder;
417  std::vector<double> _glAbscissas;
418  std::vector<double> _glWeights;
419 
420  //==================================== Variables regarding the coefficient rotation
421  bool _structureCoeffsRotated;
422 
423  //==================================== Internal functions
424  double gl20IntRR ( std::vector<double>* vals );
425  std::array<double,2> gl20IntCR ( std::vector< std::array<double,2> >* vals );
426  std::vector<double> getSingularValues ( unsigned int band, unsigned int dim, ProSHADE::ProSHADE_settings* settings );
427  std::vector< std::vector< std::vector<double> > > getSingularValuesUandVMatrix ( std::vector< std::vector<double> > mat, unsigned int dim, ProSHADE::ProSHADE_settings* settings );
428  double maxPeakHeightFromAngleAxis ( double X, double Y, double Z, double Angle, ProSHADE::ProSHADE_settings* settings );
429 
430  public:
431  //==================================== Public functions
432  ProSHADE_comparePairwise ( ProSHADE_data *cmpObj1, ProSHADE_data *cmpObj2, double mPower, std::vector<int> ignoreL, unsigned int order, ProSHADE::ProSHADE_settings* settings );
434 
435  void precomputeTrSigmaDescriptor ( );
436  void getSO3InverseMap ( ProSHADE::ProSHADE_settings* settings );
437  std::array<double,3> getEulerAngles ( ProSHADE::ProSHADE_settings* settings, double* correlation = nullptr );
438  void setEulerAngles ( double alpha, double beta, double gamma );
439  std::vector< std::array<double,8> > getSO3Peaks ( ProSHADE::ProSHADE_settings* settings, double noIQRs = 1.5, bool freeMem = true, int peakSize = 1, double realDist = 0.4, int verbose = 1 );
440  bool checkPeakExistence ( double alpha, double beta, double gamma, ProSHADE::ProSHADE_settings* settings, double passThreshold = 0.3 );
441  std::vector< std::array<double,5> > findCnSymmetry ( std::vector< std::array<double,8> > peaks,
442  ProSHADE::ProSHADE_settings* settings,
443  double axisErrorTolerance = 0.0,
444  bool freeMem = true,
445  double percAllowedToMiss = 0.33,
446  int verbose = 1 );
447  std::vector< std::array<double,5> > findCnSymmetryClear ( std::vector< std::array<double,5> > CnSymm, ProSHADE::ProSHADE_settings* settings, double maxGap = 0.2, bool *pf = nullptr );
448  std::vector< std::vector< std::array<double,6> > > findDnSymmetry ( std::vector< std::array<double,5> > CnSymm, double axisErrorTolerance = 0.1 );
449  std::vector< std::vector< std::array<double,6> > > findDnSymmetryClear ( std::vector< std::vector< std::array<double,6> > > DnSymm, ProSHADE::ProSHADE_settings* settings, double maxGap = 0.2, bool *pf = nullptr );
450  std::vector< std::vector< std::array<double,5> > > findTetrSymmetry ( std::vector< std::array<double,5> > CnSymm, double *tetrSymmPeakAvg, double axisErrorTolerance = 0.1 );
451  std::vector< std::vector< std::array<double,5> > > findOctaSymmetry ( std::vector< std::array<double,5> > CnSymm, double *octaSymmPeakAvg, double axisErrorTolerance = 0.1 );
452  std::vector< std::vector< std::array<double,5> > > findIcosSymmetry ( std::vector< std::array<double,5> > CnSymm, double *icosSymmPeakAvg, double axisErrorTolerance = 0.1 );
453  void generateWignerMatrices ( ProSHADE::ProSHADE_settings* settings );
454  double getRotCoeffDistance ( ProSHADE::ProSHADE_settings* settings );
455  double maxAvgPeakForSymmetry ( double X, double Y, double Z, double Angle, ProSHADE::ProSHADE_settings* settings );
456  void freeInvMap ( );
457  double getPeakThreshold ( );
458  void rotateStructure ( ProSHADE_data *cmpObj1, ProSHADE::ProSHADE_settings* settings,
459  std::string saveName, int verbose = 0, std::string axOrd = "xyz",
460  bool internalUse = false );
461  std::array<double,4> getTranslationFunctionMap ( ProSHADE_data *obj1, ProSHADE_data *obj2, double *ob2XMov = nullptr,
462  double *ob2YMov = nullptr, double *ob2ZMov = nullptr );
463  };
464 
473  {
474  //==================================== Friends
475  friend class ProSHADE_distances;
476 
477  private:
478  //==================================== Just save the input data and make the internal functions do the checks themselves
479  ProSHADE_data *one;
480  ProSHADE_data *two;
481  std::vector<ProSHADE_data*> *all;
482 
483  //==================================== These needs to be the same
484  std::vector<int> _lsToIgnore;
485  double _matrixPowerWeight;
486 
487  //==================================== Simple checks
488  bool _energyLevelsComputed;
489  bool _trSigmaPreComputed;
490  bool _trSigmaComputed;
491  bool _keepOrRemove;
492  bool _so3InvMapComputed;
493  bool _eulerAnglesFound;
494  bool _wignerMatricesComputed;
495  bool _fullDistComputed;
496 
497  //==================================== Result holders
498  std::vector<double> _distancesEnergyLevels;
499  std::vector<double> _distancesTraceSigma;
500  std::vector<double> _distancesFullRotation;
501 
502  //==================================== Internal value holders
503  std::vector< std::vector< std::vector< std::vector< std::array<double,2> > > > > _EMatrices;
504  std::vector<fftw_complex*> _so3InverseCoeffs;
505  std::vector< std::array<double,3> > _eulerAngles;
506  std::vector< std::vector< std::vector< std::vector< std::array<double,2> > > > > _wignerMatrices;
507  std::vector<unsigned int> _enLevelsDoNotFollow;
508  std::vector<unsigned int> _trSigmaDoNotFollow;
509 
510  //==================================== Internal functions
511  double glIntRR ( std::vector<double>* vals, double shellSep, std::vector<double>* glAbscissas, std::vector<double>* glWeights );
512  std::array<double,2> glIntCR ( std::vector< std::array<double,2> >* vals, double shellSep, std::vector<double> *glAbscissas, std::vector<double> *glWeights );
513  std::vector<double> getSingularValues ( unsigned int strNo, unsigned int band, unsigned int dim, ProSHADE::ProSHADE_settings* settings );
514 
515  public:
516  //==================================== Public functions
517  ProSHADE_compareOneAgainstAll ( ProSHADE_data *oneStr, std::vector<ProSHADE_data*> *allStrs,
518  std::vector<int> ignoreL, double matrixPowerWeight, int verbose );
519  ProSHADE_compareOneAgainstAll ( ProSHADE_data *oneStr, ProSHADE_data *twoStr, std::vector<ProSHADE_data*> *allStrs,
520  std::vector<int> ignoreL, double matrixPowerWeight, int verbose );
522 
523  std::vector<double> getEnergyLevelsDistance ( int verbose, ProSHADE::ProSHADE_settings* settings );
524  void precomputeTrSigmaDescriptor ( double shellSpacing, std::vector<unsigned int>* glIntegOrderVec, ProSHADE::ProSHADE_settings* settings );
525  std::vector<double> getTrSigmaDistance ( int verbose, ProSHADE::ProSHADE_settings* settings );
526  void getSO3InverseMap ( ProSHADE::ProSHADE_settings* settings );
527  std::vector< std::array<double,3> > getEulerAngles ( ProSHADE::ProSHADE_settings* settings );
528  void generateWignerMatrices ( ProSHADE::ProSHADE_settings* settings );
529  std::vector<double> getRotCoeffDistance ( int verbose, ProSHADE::ProSHADE_settings* settings );
530  void alignDensities ( ProSHADE::ProSHADE_settings* settings );
531  };
532 
541  {
542  private:
543  //==================================== Settings regarding resolutions
544  double mapResolution;
545  unsigned int bandwidth;
546  unsigned int glIntegOrder;
547  unsigned int theta;
548  unsigned int phi;
549  std::vector<unsigned int> bandwidthVec;
550  std::vector<unsigned int> glIntegOrderVec;
551  std::vector<unsigned int> thetaVec;
552  std::vector<unsigned int> phiVec;
553 
554  //==================================== Settings regarding B factors
555  double bFactorValue;
556  double bFactorChange;
557 
558  //==================================== Setting regarding maps and removing noise
559  double noIQRsFromMap;
560 
561  //==================================== Settings regarding concentric shells
562  double shellSpacing;
563  unsigned int manualShells;
564 
565  //==================================== Settings regarding map with phases
566  bool useCOM;
567  bool firstLineCOM;
568 
569  //==================================== Settings regarding space around the structure in lattice
570  double extraSpace;
571  std::vector<double> extraSpaceVec;
572 
573  //==================================== Settings regarding weighting the distances
574  double alpha;
575  double mPower;
576 
577  //==================================== Settings regarding bands to ignore
578  std::vector<int> ignoreLs;
579 
580  //==================================== Settings regarding the structures to use
581  std::vector <std::string> structFiles;
582 
583  //==================================== Settings regarding which distances to compute
584  bool energyLevelDist;
585  bool traceSigmaDist;
586  bool fullRotFnDist;
587 
588  //==================================== Settings regarding dealing with phases
589  bool usePhase;
590  bool saveWithAndWithout;
591 
592  //==================================== Distances holder variables
593  std::vector<double> energyLevelsDistances;
594  std::vector<double> traceSigmaDistances;
595  std::vector<double> fullRotationDistances;
596  std::vector< std::vector<double> > fragEnergyLevelsDistances;
597  std::vector< std::vector<double> > fragTraceSigmaDistances;
598  std::vector< std::vector<double> > fragfullRotationDistances;
599 
600  //==================================== Settings regarding thresholds for hierarchical distance computation
601  double enLevelsThreshold;
602  double trSigmaThreshold;
603 
604  //==================================== Pointer to the comparison object
606 
607  public:
608  //==================================== Public functions
610  ~ProSHADE_distances ( );
611  void saveDatabase ( ProSHADE::ProSHADE_settings* settings );
612  void compareAgainstDatabase ( ProSHADE::ProSHADE_settings* settings,
613  std::vector<std::string>* matchedStrNames );
614  void compareFragAgainstDatabase ( ProSHADE::ProSHADE_settings* settings,
615  std::vector<std::string>* matchedStrNames );
616  std::vector<double> getEnergyLevelsDistances ( );
617  std::vector<double> getTraceSigmaDistances ( );
618  std::vector<double> getFullRotationDistances ( );
619  std::vector< std::vector<double> > getFragEnergyLevelsDistances( );
620  std::vector< std::vector<double> > getFragTraceSigmaDistances ( );
621  std::vector< std::vector<double> > getFragFullRotationDistances( );
622  };
623 
632  {
633  private:
634  //==================================== Private variables
636  double fragBoxSize;
637  std::string mapFragName;
638  double mapFragBoxFraction;
639  double maxMapU;
640  double maxMapV;
641  double maxMapW;
642  double xRange;
643  double yRange;
644  double zRange;
645 
646  double minDensPreNorm;
647  double maxDensPreNorm;
648  double minDensPostNorm;
649  double maxDensPostNorm;
650  double postNormMean;
651  double postNormSdev;
652  double maskVolume;
653  double totalVolume;
654  double maskMean;
655  double maskSdev;
656  double maskMax;
657  double maskMin;
658  double maskDensityRMS;
659  double allDensityRMS;
660  std::array<double,3> origRanges;
661  std::array<double,3> origDims;
662  std::vector<std::string> fragFiles;
663  bool fragFilesExist;
664 
665  public:
666  //==================================== Public functions
669 
670  void printInfo ( int verbose );
671  void fragmentMap ( std::string axOrder, int verbose, ProSHADE::ProSHADE_settings* settings );
672  std::vector<std::string> getFragmentsList ( void );
673  void dealWithHalfMaps ( ProSHADE::ProSHADE_settings* settings, double *xTranslate, double *yTranslate, double *zTranslate );
674  };
675 
676 }
677 
678 
679 //============================================ Include once
680 #endif
This class deals with reading in the data and computing structure specific information including the ...
This class is responsible for reading in map and computing all its features and other possible map ma...
std::vector< std::vector< std::array< double, 6 > > > dnSymmClear
This variable holds the gap corrected Dihedral symmetry results.
unsigned int checkFileType(std::string fileName)
This function checks the input file for being either PDB or MAP formatted.
std::vector< std::array< double, 5 > > cnSymm
This variable holds the complete Cyclic symmetry results.
void printResultsClear(int verbose)
This function prints the cleared results to the screen.
std::vector< std::array< double, 5 > > generateTetrElements(std::vector< std::array< double, 5 > > symmAxes, ProSHADE::ProSHADE_settings *settings, int verbose=0)
This function generates the 12 unique tetrahedral symmetry group elements from the already detected a...
std::vector< std::array< double, 5 > > generateTetrAxes(ProSHADE_comparePairwise *cmpObj, ProSHADE::ProSHADE_settings *settings, std::vector< std::array< double, 5 > > tetrSymm, std::vector< std::array< double, 5 > > allCs, double axisErrorTolerance=0.1, int verbose=0)
This function generates all the tetrahedral symmetry axes from a pair detected by findTetrSymmetry fu...
This namespace contains all the internal objects and their forward declarations.
std::vector< std::array< double, 5 > > icosAxes
This variable holds the 31 unique axes of the octahedral symmetry, if they are found.
std::vector< std::array< double, 5 > > generateOctaElements(std::vector< std::array< double, 5 > > symmAxes, ProSHADE::ProSHADE_settings *settings, int verbose=0)
This function generates the 24 octahedral symmetry group elements from the axes generated by generate...
std::vector< std::array< double, 5 > > generateOctaAxes(ProSHADE_comparePairwise *cmpObj, ProSHADE::ProSHADE_settings *settings, std::vector< std::array< double, 5 > > octaSymm, std::vector< std::array< double, 5 > > allCs, double axisErrorTolerance=0.1, int verbose=0)
This function generates all the ctahedral symmetry elements from a pair detected by findOctaSymmetry ...
This is the executive class responsible for comparing strictly two structures.
std::vector< std::array< double, 5 > > generateIcosAxes(ProSHADE_comparePairwise *cmpObj, ProSHADE::ProSHADE_settings *settings, std::vector< std::array< double, 5 > > icosSymm, std::vector< std::array< double, 5 > > allCs, double axisErrorTolerance=0.1, int verbose=0)
This function generates all the icosahedral symmetry elements from a pair detected by findIcosSymmetr...
std::vector< std::array< double, 5 > > cnSymmClear
This variable holds the gap corrected Cyclic symmetry results.
std::vector< std::array< double, 5 > > generateIcosElements(std::vector< std::array< double, 5 > > symmAxes, ProSHADE::ProSHADE_settings *settings, int verbose=0)
This function generates the 60 icosahedral symmetry group elements from the axes generated by generat...
void printResultsClearHTML(ProSHADE::ProSHADE_settings *settings)
This function prints the cleared results to the HTML file.
bool inputStructureDataType
This variable tells whether input data type is PDB or not.
std::vector< std::array< double, 5 > > getCSymmetries()
This function gives the user programmatical access to the detected C-symmetries.
std::vector< std::array< double, 8 > > getRotFnPeaks()
This function gives the user programmatical access to the symmetry peaks of the symmetry detection pr...
void printResultsRequestHTML(std::string symmetryType, unsigned int symmetryFold, ProSHADE::ProSHADE_settings *settings)
This function prints the cleared results to the HTML report file.
This is the class which computes the symmetry in a single structure.
ProSHADE_symmetry()
Contructor for the ProSHADE_symmetry class.
This is the executive class responsible for comparing two or more structures.
std::vector< std::array< double, 5 > > tetrAxes
This variable holds the 7 unique axes of the tetrahedral symmetry, if they are found.
This is the executive class for computing distances between two or more structures.
This class stores all the settings and is passed to the executive classes instead of multitude of par...
Definition: ProSHADE.h:74
std::vector< std::array< double, 5 > > octaAxes
This variable holds the 13 unique axes of the octahedral symmetry, if they are found.
std::vector< std::vector< std::array< double, 6 > > > getDSymmetries()
This function gives the user programmatical access to the detected D-symmetries.
void printResultsRequest(std::string symmetryType, unsigned int symmetryFold, int verbose)
This function prints the cleared results to the screen.
std::vector< std::vector< std::array< double, 6 > > > dnSymm
This variable holds the complete Dihedral symmetry results.