ProSHADE  0.6.5 (NOV 2018)
Protein Shape Descriptors and Symmetry Detection
All Classes Namespaces Files Functions Variables Pages
getDistances.cpp
Go to the documentation of this file.
1 
15 //============================================ ProSHADE
16 #include "../proshade/ProSHADE.h"
17 
18 //============================================ Main
19 int main ( int argc,
20  char **argv )
21 {
22  //======================================== Default settings
24 
25  // ... Settings regarding resolutions
26  setUp->mapResolution = 8.0;
27  setUp->bandwidth = 0;
28  setUp->glIntegOrder = 0;
29  setUp->theta = 0;
30  setUp->phi = 0;
31 
32  // ... Settings regarding B factors
33  setUp->bFactorValue = 80.0;
34  setUp->bFactorChange = 0.0;
35 
36  // ... Settings regarding concentric shells
37  setUp->shellSpacing = 0.0;
38  setUp->manualShells = 0;
39 
40  // ... Settings regarding phase
41  setUp->usePhase = true;
42 
43  // ... Settings regarding map with phases
44  setUp->useCOM = true;
45  setUp->maskBlurFactor = 500.0;
46  setUp->maskBlurFactorGiven = false;
47 
48  // ... Settings regarding weighting the distances
49  setUp->alpha = 1.0;
50  setUp->mPower = 1.0;
51 
52  // ... Settings regarding bands to be ignored
53  std::vector<int> lsToIgnore;
54  lsToIgnore.emplace_back ( 0 );
55  setUp->ignoreLs = lsToIgnore;
56 
57  // ... Settings regarding which distances to compute
58  setUp->energyLevelDist = true;
59  setUp->traceSigmaDist = true;
60  setUp->fullRotFnDist = true;
61 
62  // ... Settings regarding hierarchical distance computation
63  setUp->enLevelsThreshold = -999.9;
64  setUp->trSigmaThreshold = -999.9;
65 
66  // ... Settings regarding the task
67  setUp->taskToPerform = ProSHADE::Distances;
68 
69  // ... Settings regarding loudness
70  setUp->verbose = -1;
71  setUp->htmlReport = false;
72 
73  //======================================== Get files
74  if ( argc < 3 )
75  {
76  std::cout << std::endl << "Usage: getDistances [filename1] [filename2] ... [filenameX] to get distances from [filename1] to all other files. Minimum of two files." << std::endl << std::endl;
77  exit ( 0 );
78  }
79  else
80  {
81  for ( unsigned int iter = 1; iter < static_cast<unsigned int> ( argc ); iter++ )
82  {
83  setUp->structFiles.emplace_back ( std::string ( argv[iter] ) );
84  }
85  }
86 
87  //======================================== Run ProSHADE
88  ProSHADE::ProSHADE *run = new ProSHADE::ProSHADE ( setUp );
89 
90  //======================================== Get results
91  std::vector<double> crossCorrDists = run->getCrossCorrDists ( );
92  std::vector<double> traceSigmaDists = run->getTraceSigmaDists ( );
93  std::vector<double> rotFunDists = run->getRotFunctionDists ( );
94 
95  //======================================== Print results
96  printf ( "ProSHADE library version: %s\n\n", run->getProSHADEVersion().c_str() );
97 
98  if ( crossCorrDists.size() > 0 )
99  {
100  printf ( "Energy Level Descriptor distances : %+.5f", crossCorrDists.at(0) );
101  for ( unsigned int iter = 1; iter < static_cast<unsigned int> ( crossCorrDists.size() ); iter++ )
102  {
103  printf ( "\t%+.4f", crossCorrDists.at(iter) );
104  }
105  printf ( "\n" );
106  }
107  else
108  {
109  std::cerr << "!!! Error !!! ProSHADE failed to obtain distances, please see the standard error stream for details. Did you set the energyLevelDist settings to false? Terminating..." << std::endl;
110  exit ( -1 );
111  }
112 
113  if ( traceSigmaDists.size() > 0 )
114  {
115  printf ( "Trace Sigma Descriptor distances : %+.5f", traceSigmaDists.at(0) );
116  for ( unsigned int iter = 1; iter < static_cast<unsigned int> ( traceSigmaDists.size() ); iter++ )
117  {
118  printf ( "\t%+.4f", traceSigmaDists.at(iter) );
119  }
120  printf ( "\n" );
121  }
122  else
123  {
124  std::cerr << "!!! Error !!! ProSHADE failed to obtain distances, please see the standard error stream for details. Did you set the traceSigmaDist settings to false? Terminating..." << std::endl;
125  exit ( -1 );
126  }
127 
128  if ( rotFunDists.size() > 0 )
129  {
130  printf ( "Rotation Function Descriptor distances : %+.5f", rotFunDists.at(0) );
131  for ( unsigned int iter = 1; iter < static_cast<unsigned int> ( rotFunDists.size() ); iter++ )
132  {
133  printf ( "\t%+.4f", rotFunDists.at(iter) );
134  }
135  printf ( "\n" );
136  }
137  else
138  {
139  std::cerr << "!!! Error !!! ProSHADE failed to obtain distances, please see the standard error stream for details. Did you set the fullRotFnDist settings to false? Terminating..." << std::endl;
140  exit ( -1 );
141  }
142 
143  for ( unsigned int iter = 1; iter < static_cast<unsigned int> ( setUp->structFiles.size() ); iter++ )
144  {
145  printf ( "Matching structure names : %40s vs. %-40s\n", setUp->structFiles.at(0).c_str(), setUp->structFiles.at(iter).c_str() );
146  }
147  std::cout << std::endl;
148 
149  //======================================== Free memory
150  delete setUp;
151  delete run;
152 
153  //======================================== Done
154  return 0;
155 }
156 
double mapResolution
This is the internal resolution at which the calculations are done, not necessarily the resolution of...
Definition: ProSHADE.h:78
std::vector< double > getRotFunctionDists(void)
Accessor function for the rotation function based distances vector.
Definition: ProSHADE.cpp:3011
unsigned int theta
This parameter is the longitude of the spherical grid mapping. It should be 2 * bandwidth unless ther...
Definition: ProSHADE.h:83
double bFactorValue
This is the value to which all B-factors of PDB files will be changed to.
Definition: ProSHADE.h:88
bool usePhase
Here the user can decide whether to use phase information or whether to ignore it completely...
Definition: ProSHADE.h:101
unsigned int bandwidth
This parameter determines the angular resolution of the spherical harmonics decomposition.
Definition: ProSHADE.h:80
bool htmlReport
Should HTML report for the run be created?
Definition: ProSHADE.h:186
int verbose
Should the software report on the progress, or just be quiet? Value between 0 (quiet) and 4 (loud) ...
Definition: ProSHADE.h:189
std::string getProSHADEVersion(void)
Miscellanous function allowing the user to get the ProSHADE version.
Definition: ProSHADE.cpp:2910
double trSigmaThreshold
All structure pairs with trace sigma descriptor value less than this will not be subjected to any fur...
Definition: ProSHADE.h:138
double alpha
This parameter determines the power to which the |F|&#39;s should be raised.
Definition: ProSHADE.h:113
bool fullRotFnDist
Should the full rotation function distances descriptor be computed.
Definition: ProSHADE.h:134
bool maskBlurFactorGiven
Was a specific value of the blurring factor requested by the user?
Definition: ProSHADE.h:148
std::vector< int > ignoreLs
This vector lists all the bandwidth values which should be ignored and not part of the computations...
Definition: ProSHADE.h:117
double shellSpacing
This parameter determines how far the radial shells should be from each other.
Definition: ProSHADE.h:96
std::vector< double > getTraceSigmaDists(void)
Accessor function for the trace sigma distances vector.
Definition: ProSHADE.cpp:2985
bool traceSigmaDist
Should the trace sigma distances descriptor be computed.
Definition: ProSHADE.h:133
unsigned int phi
This parameter is the latitudd of the spherical grid mapping. It should be 2 * bandwidth unless there...
Definition: ProSHADE.h:84
bool energyLevelDist
Should the energy level distances descriptor be computed.
Definition: ProSHADE.h:132
double bFactorChange
This value will be used to change the B-factors if required by the user.
Definition: ProSHADE.h:89
This class stores all the settings and is passed to the executive classes instead of multitude of par...
Definition: ProSHADE.h:74
double mPower
This parameter determines the scaling for trace sigma descriptor.
Definition: ProSHADE.h:114
unsigned int manualShells
Should the user require so, the maximum number of radial shells can be set.
Definition: ProSHADE.h:98
std::vector< double > getCrossCorrDists(void)
Accessor function for the cross-correlation distances vector.
Definition: ProSHADE.cpp:2959
bool useCOM
Should the Centre of Mass (COM) be used to center the structure in the cell?
Definition: ProSHADE.h:105
double maskBlurFactor
The is the amount of blurring to be used to create masks for maps.
Definition: ProSHADE.h:147
std::vector< std::string > structFiles
This vector should contain all the structures that are being dealt with, but this does not yet work! ...
Definition: ProSHADE.h:120
double enLevelsThreshold
All structure pairs with energy level descriptor value less than this will not be subjected to any fu...
Definition: ProSHADE.h:137
Task taskToPerform
This custom type variable determines which task to perfom (i.e. symmetry detection, distances computation or map features extraction).
Definition: ProSHADE.h:141
unsigned int glIntegOrder
This parameter controls the Gauss-Legendre integration order and so the radial resolution.
Definition: ProSHADE.h:82