ProSHADE  0.6.5 (NOV 2018)
Protein Shape Descriptors and Symmetry Detection
All Classes Namespaces Files Functions Variables Pages
buildDatabase.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 = 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  // ... 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 phase
44  setUp->usePhase = true;
45 
46  // ... Settings regarding map with phases
47  setUp->useCOM = true;
48  setUp->firstLineCOM = false;
49 
50  // ... Settings regarding space around the structure in lattice
51  setUp->extraSpace = 8.0;
52 
53  // ... Settings regarding bands to be ignored
54  std::vector<int> lsToIgnore;
55  lsToIgnore.emplace_back ( 0 );
56  setUp->ignoreLs = lsToIgnore;
57 
58  // ... Settings regarding the task
59  setUp->taskToPerform = ProSHADE::BuildDB;
60 
61  // ... Settings regarding where and if to save the clear map
62  setUp->clearMapData = true;
63 
64  // ... Settings regarding loudness
65  setUp->verbose = -1;
66  setUp->htmlReport = false;
67 
68  // ... Settings regarding the database
69  setUp->databaseName = "";
70 
71  //======================================== Get files
72  if ( argc < 3 )
73  {
74  std::cout << std::endl << "Usage: buildDatabase [dbName] [filename1] [filename2] ... [filenameX] build the database in [dbName] using files [filename1] to [filenameX]. Minimum of two files." << std::endl << std::endl;
75  exit ( 0 );
76  }
77  else
78  {
79  setUp->databaseName = std::string ( argv[1] ) ;
80  for ( unsigned int iter = 2; iter < static_cast<unsigned int> ( argc ); iter++ )
81  {
82  setUp->structFiles.emplace_back ( std::string ( argv[iter] ) );
83  }
84  }
85 
86  //======================================== Run ProSHADE
87  ProSHADE::ProSHADE *run = new ProSHADE::ProSHADE ( setUp );
88 
89  //======================================== Print results
90  printf ( "ProSHADE library version: %s\n\n", run->getProSHADEVersion().c_str() );
91 
92  //======================================== Free memory
93  delete setUp;
94  delete run;
95 
96  //======================================== Done
97  return 0;
98 }
99 
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
bool clearMapData
This value is used to decide whether the input maps should be cleared again, or not.
Definition: ProSHADE.h:146
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::string databaseName
The name of the bin file to which the database should be saved.
Definition: ProSHADE.h:156
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
bool firstLineCOM
This is a special option for metal detection, please leave false.
Definition: ProSHADE.h:106
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
unsigned int phi
This parameter is the latitudd of the spherical grid mapping. It should be 2 * bandwidth unless there...
Definition: ProSHADE.h:84
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
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
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
unsigned int glIntegOrder
This parameter controls the Gauss-Legendre integration order and so the radial resolution.
Definition: ProSHADE.h:82