ProSHADE  0.6.5 (NOV 2018)
Protein Shape Descriptors and Symmetry Detection
All Classes Namespaces Files Functions Variables Pages
getSymmetry.cpp
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 = 6.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  // ... Setting regarding maps and removing noise
37  setUp->noIQRsFromMap = 4.0;
38 
39  // ... Settings regarding concentric shells
40  setUp->shellSpacing = 0.0;
41  setUp->manualShells = 0;
42 
43  // ... Settings regarding map with phases
44  setUp->useCOM = false;
45  setUp->maskBlurFactor = 500.0;
46  setUp->maskBlurFactorGiven = false;
47 
48  // ... Settings regarding space around the structure in lattice
49  setUp->extraSpace = 8.0;
50 
51  // ... Settings regarding peak detection
52  setUp->peakHeightNoIQRs = 3.0;
53  setUp->peakDistanceForReal = 0.20;
54  setUp->peakSurroundingPoints = 1;
55 
56  // ... Settings regarding tolerances
57  setUp->aaErrorTolerance = 0.1;
58  setUp->symGapTolerance = 0.3;
59 
60  // ... Settings regarding the task
61  setUp->taskToPerform = ProSHADE::Symmetry;
62 
63  // ... Settings regarding the symmetry type required
64  setUp->symmetryFold = 0;
65  setUp->symmetryType = "";
66 
67  // ... Settings regarding loudness
68  setUp->verbose = -1;
69  setUp->htmlReport = false;
70 
71  //======================================== Get files
72  if ( argc != 2 )
73  {
74  std::cout << std::endl << "Usage: getSymmetry [filename1] to get symmetry information for structure in [filename1]." << std::endl << std::endl;
75  exit ( 0 );
76  }
77  else
78  {
79  setUp->structFiles.emplace_back ( std::string ( argv[1] ) );
80  }
81 
82  //======================================== Run ProSHADE
83  ProSHADE::ProSHADE *run = new ProSHADE::ProSHADE ( setUp );
84 
85  //======================================== Get results
86  std::vector< std::array<double,5> > cSyms = run->getCyclicSymmetries ( );
87  std::vector< std::vector< std::array<double,6> > > dSyms = run->getDihedralSymmetries ( );
88  std::vector< std::array<double,5> > tetSym= run->getTetrahedralSymmetries ( );
89  std::vector< std::array<double,5> > octSym= run->getOctahedralSymmetries ( );
90  std::vector< std::array<double,5> > icoSym= run->getIcosahedralSymmetries ( );
91  std::vector< std::array<double,5> > symRecom = run->getRecommendedSymmetry ( );
92  std::vector< std::array<double,5> > symSpecElems = run->getSpecificSymmetryElements ( "C", 4 );
93 
94  //======================================== Print results
95  printf ( "ProSHADE library version: %s\n\n", run->getProSHADEVersion().c_str() );
96 
97  printf ( "Cyclic symmetry axes detected:\n" );
98  printf ( "-----------------------------------------------------------\n" );
99  printf ( "Symmetry Fold x y z Angle Peak\n" );
100  printf ( " Type height\n" );
101  printf ( "-----------------------------------------------------------\n" );
102  for ( unsigned int iter = 0; iter < static_cast<unsigned int> ( cSyms.size() ); iter++ )
103  {
104  printf ( " C %d %+.2f %+.2f %+.2f 2pi / %d %+.3f\n", static_cast<int> ( cSyms.at(iter)[0] ), cSyms.at(iter)[1], cSyms.at(iter)[2], cSyms.at(iter)[3], static_cast<int> ( cSyms.at(iter)[0] ), static_cast<double> ( cSyms.at(iter)[4] ) );
105  }
106  printf ( "\n" );
107 
108  printf ( "Dihedral symmetry axes table:\n" );
109  printf ( "-----------------------------------------------------------\n" );
110  printf ( "Symmetry Fold x y z Angle Peak\n" );
111  printf ( " Type height\n" );
112  printf ( "-----------------------------------------------------------\n" );
113  for ( unsigned int iter = 0; iter < static_cast<unsigned int> ( dSyms.size() ); iter++ )
114  {
115  printf ( " D %d %+.2f %+.2f %+.2f 2pi / %d %+.3f\n", static_cast<int> ( dSyms.at(iter).at(0)[0] ), dSyms.at(iter).at(0)[1], dSyms.at(iter).at(0)[2], dSyms.at(iter).at(0)[3], static_cast<int> ( dSyms.at(iter).at(0)[0] ), dSyms.at(iter).at(0)[4] );
116  for ( unsigned int it = 1; it < static_cast<unsigned int> ( dSyms.at(iter).size() ); it++ )
117  {
118  printf ( " %d %+.2f %+.2f %+.2f 2pi / %d %+.3f\n", static_cast<int> ( dSyms.at(iter).at(it)[0] ), dSyms.at(iter).at(it)[1], dSyms.at(iter).at(it)[2], dSyms.at(iter).at(it)[3], static_cast<int> ( dSyms.at(iter).at(it)[0] ), dSyms.at(iter).at(it)[4] );
119  }
120 
121  }
122  printf ( "\n" );
123 
124  printf ( "Tetrahedral symmetry axes table:\n" );
125  printf ( "-----------------------------------------------------------\n" );
126  printf ( "Symmetry Fold x y z Angle Peak\n" );
127  printf ( " Type height\n" );
128  printf ( "-----------------------------------------------------------\n" );
129  if ( tetSym.size() > 0 )
130  {
131  printf ( " T %d %+.2f %+.2f %+.2f 2pi / %d %+.3f\n", static_cast<int> ( tetSym.at(0)[0] ), tetSym.at(0)[1], tetSym.at(0)[2], tetSym.at(0)[3], static_cast<int> ( tetSym.at(0)[0] ), static_cast<double> ( tetSym.at(0)[4] ) );
132  for ( unsigned int iter = 1; iter < static_cast<unsigned int> ( tetSym.size() ); iter++ )
133  {
134  printf ( " %d %+.2f %+.2f %+.2f 2pi / %d %+.3f\n", static_cast<int> ( tetSym.at(iter)[0] ), tetSym.at(iter)[1], tetSym.at(iter)[2], tetSym.at(iter)[3], static_cast<int> ( tetSym.at(iter)[0] ), static_cast<double> ( tetSym.at(iter)[4] ) );
135  }
136  }
137  printf ( "\n" );
138 
139  printf ( "Octahedral symmetry axes table:\n" );
140  printf ( "-----------------------------------------------------------\n" );
141  printf ( "Symmetry Fold x y z Angle Peak\n" );
142  printf ( " Type height\n" );
143  printf ( "-----------------------------------------------------------\n" );
144  if ( octSym.size() > 0 )
145  {
146  printf ( " O %d %+.2f %+.2f %+.2f 2pi / %d %+.3f\n", static_cast<int> ( octSym.at(0)[0] ), octSym.at(0)[1], octSym.at(0)[2], octSym.at(0)[3], static_cast<int> ( octSym.at(0)[0] ), static_cast<double> ( octSym.at(0)[4] ) );
147  for ( unsigned int iter = 1; iter < static_cast<unsigned int> ( octSym.size() ); iter++ )
148  {
149  printf ( " %d %+.2f %+.2f %+.2f 2pi / %d %+.3f\n", static_cast<int> ( octSym.at(iter)[0] ), octSym.at(iter)[1], octSym.at(iter)[2], octSym.at(iter)[3], static_cast<int> ( octSym.at(iter)[0] ), static_cast<double> ( octSym.at(iter)[4] ) );
150  }
151  }
152  printf ( "\n" );
153 
154  printf ( "Icosahedral symmetry axes table:\n" );
155  printf ( "-----------------------------------------------------------\n" );
156  printf ( "Symmetry Fold x y z Angle Peak\n" );
157  printf ( " Type height\n" );
158  printf ( "-----------------------------------------------------------\n" );
159  if ( icoSym.size() > 0 )
160  {
161  printf ( " I %d %+.2f %+.2f %+.2f 2pi / %d %+.3f\n", static_cast<int> ( icoSym.at(0)[0] ), icoSym.at(0)[1], icoSym.at(0)[2], icoSym.at(0)[3], static_cast<int> ( icoSym.at(0)[0] ), static_cast<double> ( icoSym.at(0)[4] ) );
162  for ( unsigned int iter = 1; iter < static_cast<unsigned int> ( icoSym.size() ); iter++ )
163  {
164  printf ( " %d %+.2f %+.2f %+.2f 2pi / %d %+.3f\n", static_cast<int> ( icoSym.at(iter)[0] ), icoSym.at(iter)[1], icoSym.at(iter)[2], icoSym.at(iter)[3], static_cast<int> ( icoSym.at(iter)[0] ), static_cast<double> ( icoSym.at(iter)[4] ) );
165  }
166  }
167  printf ( "\n" );
168 
169  printf ( "-----------------------------------------------------------\n" );
170  printf ( "-----------------------------------------------------------\n" );
171  printf ( "Recommended symmetry: " );
172  if ( icoSym.size() > 0 ) { printf ( "ICOSAHEDRAL\n" ); }
173  else if ( octSym.size() > 0 ) { printf ( "OCTAHEDRAL\n" ); }
174  else if ( tetSym.size() > 0 ) { printf ( "TETRAHEDRAL\n" ); }
175  else if ( dSyms.size() > 0 ) { printf ( "DIHEDRAL\n" ); }
176  else if ( cSyms.size() > 0 ) { printf ( "CYCLIC\n" ); }
177  printf ( "-----------------------------------------------------------\n" );
178  printf ( "-----------------------------------------------------------\n\n" );
179  printf ( "Recommended symmetry elements table:\n" );
180  printf ( "-----------------------------------------------------------\n" );
181  printf ( "Symmetry x y z Angle \n" );
182  printf ( " Type (rad) \n" );
183  printf ( "-----------------------------------------------------------\n" );
184  if ( symRecom.size() > 0 )
185  {
186  printf ( " E %+.2f %+.2f %+.2f %+.3f \n", symRecom.at(0)[1], symRecom.at(0)[2], symRecom.at(0)[3], symRecom.at(0)[4] );
187  for ( unsigned int iter = 1; iter < static_cast<unsigned int> ( symRecom.size() ); iter++ )
188  {
189  printf ( " C%d %+.2f %+.2f %+.2f %+.3f \n", static_cast<int> ( symRecom.at(iter)[0] ), symRecom.at(iter)[1], symRecom.at(iter)[2], symRecom.at(iter)[3], symRecom.at(iter)[4] );
190  }
191  }
192  printf ( "\n" );
193 
194  printf ( "Requested symmetry elements table:\n" );
195  printf ( "-----------------------------------------------------------\n" );
196  printf ( "Symmetry x y z Angle \n" );
197  printf ( " Type (rad) \n" );
198  printf ( "-----------------------------------------------------------\n" );
199  if ( symSpecElems.size() > 0 )
200  {
201  printf ( " E %+.2f %+.2f %+.2f %+.3f \n", symSpecElems.at(0)[1], symSpecElems.at(0)[2], symSpecElems.at(0)[3], symSpecElems.at(0)[4] );
202  for ( unsigned int iter = 1; iter < static_cast<unsigned int> ( symSpecElems.size() ); iter++ )
203  {
204  printf ( " C%d %+.2f %+.2f %+.2f %+.3f \n", static_cast<int> ( symSpecElems.at(iter)[0] ), symSpecElems.at(iter)[1], symSpecElems.at(iter)[2], symSpecElems.at(iter)[3], symSpecElems.at(iter)[4] );
205  }
206  }
207 
208  //======================================== Free memory
209  delete setUp;
210  delete run;
211 
212  //======================================== Done
213  return 0;
214 }
215 
double aaErrorTolerance
The tolerance parameter on matching axes for the angle-axis representation of rotations.
Definition: ProSHADE.h:128
std::vector< std::array< double, 5 > > getRecommendedSymmetry(void)
Accessor function for the recommended symmetry elements list.
Definition: ProSHADE.cpp:3419
std::string symmetryType
The required symmetry type. If no symmetry is required, leave empty. Possible values are: C...
Definition: ProSHADE.h:163
double mapResolution
This is the internal resolution at which the calculations are done, not necessarily the resolution of...
Definition: ProSHADE.h:78
double noIQRsFromMap
This is the number of interquartile distances from mean that is used to threshold the map masking...
Definition: ProSHADE.h:93
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
std::vector< std::array< double, 5 > > getIcosahedralSymmetries(void)
Accessor function for the icosahedral symmetries list.
Definition: ProSHADE.cpp:3354
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
std::vector< std::array< double, 5 > > getCyclicSymmetries(void)
Accessor function for the cyclic symmetries list.
Definition: ProSHADE.cpp:3037
bool maskBlurFactorGiven
Was a specific value of the blurring factor requested by the user?
Definition: ProSHADE.h:148
int peakSurroundingPoints
For a peak to exist, how many points in every direction need to be smalled than the middle value...
Definition: ProSHADE.h:125
double shellSpacing
This parameter determines how far the radial shells should be from each other.
Definition: ProSHADE.h:96
std::vector< std::array< double, 5 > > getSpecificSymmetryElements(std::string symType, int symFold=0)
Accessor function for the given symmetry elements list.
Definition: ProSHADE.cpp:3635
double peakDistanceForReal
Threshold for determining &#39;missing peaks&#39; existence.
Definition: ProSHADE.h:124
unsigned int phi
This parameter is the latitudd of the spherical grid mapping. It should be 2 * bandwidth unless there...
Definition: ProSHADE.h:84
unsigned int symmetryFold
The required fold of the sought symmetry. Applicable to C and D symmetries, otherwise leave 0...
Definition: ProSHADE.h:162
std::vector< std::array< double, 5 > > getTetrahedralSymmetries(void)
Accessor function for the tetrahedral symmetries list.
Definition: ProSHADE.cpp:3226
double bFactorChange
This value will be used to change the B-factors if required by the user.
Definition: ProSHADE.h:89
double peakHeightNoIQRs
How many interquartile ranges should be used to distinguish &#39;false&#39; peaks from the true ones...
Definition: ProSHADE.h:123
std::vector< std::array< double, 5 > > getOctahedralSymmetries(void)
Accessor function for the octahedral symmetries list.
Definition: ProSHADE.cpp:3290
double symGapTolerance
For C-symmetries - if there are many, only those with average peak height - parameter * top symmetry ...
Definition: ProSHADE.h:129
This class stores all the settings and is passed to the executive classes instead of multitude of par...
Definition: ProSHADE.h:74
unsigned int manualShells
Should the user require so, the maximum number of radial shells can be set.
Definition: ProSHADE.h:98
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
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
double extraSpace
What should be the distance added on both sides to the structure, so that the next cell density would...
Definition: ProSHADE.h:109
std::vector< std::vector< std::array< double, 6 > > > getDihedralSymmetries(void)
Accessor function for the dihedral symmetries list. These will be sorted using the fold and peak heig...
Definition: ProSHADE.cpp:3173
unsigned int glIntegOrder
This parameter controls the Gauss-Legendre integration order and so the radial resolution.
Definition: ProSHADE.h:82