ProSHADE  0.6.6 (DEC 2018)
Protein Shape Descriptors and Symmetry Detection
mapSearch.cpp
Go to the documentation of this file.
1 
20 //============================================ ProSHADE
21 #include "../proshade/ProSHADE.h"
22 #include "../proshade/ProSHADE_internal.h"
23 
24 //============================================ Main
25 int main ( int argc, char **argv )
26 {
27  //======================================== Initialise variables
29  setUp->verbose = -1;
30  setUp->extraSpace = 0.0;
31  setUp->volumeTolerance = 10.0;
32  setUp->useCOM = false;
33 
34  std::string mapName = "../xx_searchMap_data/02/emd_7780.map";
35  setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainA.pdb" );
36 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainB.pdb" );
37 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainC.pdb" );
38 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainD.pdb" );
39 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainE.pdb" );
40 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainG.pdb" );
41 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainH.pdb" );
42 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainI.pdb" );
43 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainJ.pdb" );
44 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainK.pdb" );
45 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainL.pdb" );
46 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainM.pdb" );
47 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainN.pdb" );
48 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainO.pdb" );
49 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainP.pdb" );
50 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainQ.pdb" );
51 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainS.pdb" );
52 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainT.pdb" );
53 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainU.pdb" );
54 // setUp->structFiles.emplace_back ( "../xx_searchMap_data/02/ChainV.pdb" );
55 
56  std::vector<std::array<double,2> > strRadii;
57  std::vector<double> shSpacing ( setUp->structFiles.size(), setUp->shellSpacing );
58  std::vector<unsigned int> bandVec ( setUp->structFiles.size(), setUp->bandwidth );
59  std::vector<unsigned int> thetaVec ( setUp->structFiles.size(), setUp->theta );
60  std::vector<unsigned int> phiVec ( setUp->structFiles.size(), setUp->phi );
61  std::vector<unsigned int> glIntegVec ( setUp->structFiles.size(), setUp->glIntegOrder );
62  std::vector<double> extraVec ( setUp->structFiles.size(), setUp->extraSpace );
63  std::array<double,2> hlpArr;
64 
65  double maxRange = 0.0;
66  double secRange = 0.0;
67 
68  //======================================== For each input structure
69  for ( unsigned int strIt = 0; strIt < static_cast<unsigned int> ( setUp->structFiles.size ( ) ); strIt++ )
70  {
71  //==================================== Initialise structure reading
73 
74  //==================================== Read in the structures
75  int fileType = ProSHADE_internal::checkFileType ( setUp->structFiles.at ( strIt ) );
76  if ( fileType == 1 )
77  {
78  str->getDensityMapFromPDB ( setUp->structFiles.at ( strIt ),
79  &shSpacing.at ( strIt ),
80  setUp->mapResolution,
81  &bandVec.at ( strIt ),
82  &thetaVec.at ( strIt ),
83  &phiVec.at ( strIt ),
84  &glIntegVec.at ( strIt ),
85  &extraVec.at ( strIt ),
86  setUp->mapResDefault,
87  setUp,
88  setUp->bFactorValue,
89  false,
90  setUp->overlayDefaults );
91  }
92  else
93  {
94  str->getDensityMapFromMAP ( setUp->structFiles.at ( strIt ),
95  &shSpacing.at ( strIt ),
96  setUp->mapResolution,
97  &bandVec.at ( strIt ),
98  &thetaVec.at ( strIt ),
99  &phiVec.at ( strIt ),
100  &glIntegVec.at ( strIt ),
101  &extraVec.at ( strIt ),
102  setUp->mapResDefault,
103  setUp->rotChangeDefault,
104  setUp,
105  setUp->overlayDefaults );
106  }
107 
108  //==================================== Find input structure radius
109  if ( fileType == 1 )
110  {
111  hlpArr[0] = strIt;
112  hlpArr[1] = setUp->modelRadius ( setUp->structFiles.at ( strIt ) ) * 1.1;
113  strRadii.emplace_back ( hlpArr );
114  }
115  else
116  {
117  maxRange = std::max ( str->getMapXRange ( ), std::max ( str->getMapYRange ( ), str->getMapZRange ( ) ) );
118  if ( maxRange == str->getMapXRange ( ) )
119  {
120  secRange = std::max ( str->getMapYRange ( ), str->getMapZRange ( ) );
121  }
122  else if ( maxRange == str->getMapYRange ( ) )
123  {
124  secRange = std::max ( str->getMapXRange ( ), str->getMapZRange ( ) );
125  }
126  else
127  {
128  secRange = std::max ( str->getMapXRange ( ), str->getMapYRange ( ) );
129  }
130 
131  hlpArr[0] = strIt;
132  hlpArr[1] = sqrt ( pow ( maxRange, 2.0 ) + pow ( secRange, 2.0 ) );
133  strRadii.emplace_back ( hlpArr );
134  }
135 
136  //==================================== Free memory
137  delete str;
138  }
139 
140  //======================================== Sort by radii
141  std::sort ( strRadii.begin (), strRadii.end(), [](const std::array<double,2>& a, const std::array<double,2>& b) { return a[1] > b[1]; } );
142 
143  //======================================== Find settings max's and min's
144  double minSpacing = std::numeric_limits<double>::infinity ( );
145  unsigned int maxBand = 0;
146  unsigned int maxTheta = 0;
147  unsigned int maxPhi = 0;
148  unsigned int maxGlInt = 0;
149  double maxExtraSpace = 0.0;
150  for ( unsigned int iter = 0; iter < static_cast<unsigned int> ( setUp->structFiles.size ( ) ); iter++ )
151  {
152  if ( shSpacing.at ( iter ) < minSpacing ) { minSpacing = shSpacing.at ( iter ); }
153  if ( bandVec.at ( iter ) > maxBand ) { maxBand = bandVec.at ( iter ); }
154  if ( thetaVec.at ( iter ) > maxTheta ) { maxTheta = thetaVec.at ( iter ); }
155  if ( phiVec.at ( iter ) > maxPhi ) { maxPhi = phiVec.at ( iter ); }
156  if ( glIntegVec.at ( iter ) > maxGlInt ) { maxGlInt = glIntegVec.at ( iter ); }
157  if ( extraVec.at ( iter ) > maxExtraSpace) { maxExtraSpace = extraVec.at ( iter ); }
158  }
159 
160  double minSpacingCopy = minSpacing;
161  unsigned int maxBandCopy = maxBand;
162  unsigned int maxThetaCopy = maxTheta;
163  unsigned int maxPhiCopy = maxPhi;
164  unsigned int maxGlIntCopy = maxGlInt;
165  double maxExtraSpaceCopy = maxExtraSpace;
166 
167  //======================================== Build the database
168  setUp->taskToPerform = ProSHADE::BuildDB;
169  setUp->databaseName = "testDB.bin";
170 
171  setUp->shellSpacing = minSpacing;
172  setUp->bandwidth = maxBand;
173  setUp->theta = maxTheta;
174  setUp->phi = maxPhi;
175  setUp->glIntegOrder = maxGlInt;
176  setUp->extraSpace = maxExtraSpace;
177 
178  setUp->clearMapData = false;
179  setUp->saveWithAndWithout = true;
180 
181  ProSHADE::ProSHADE *run = new ProSHADE::ProSHADE ( setUp );
182  delete run;
183 
184  //======================================== Read in the map
186  inMap->getDensityMapFromMAP ( mapName,
187  &minSpacingCopy,
188  setUp->mapResolution,
189  &maxBandCopy,
190  &maxThetaCopy,
191  &maxPhiCopy,
192  &maxGlIntCopy,
193  &maxExtraSpaceCopy,
194  setUp->mapResDefault,
195  setUp->rotChangeDefault,
196  setUp,
197  setUp->overlayDefaults );
198 
199  //======================================== Map the query to the shells
200  inMap->keepPhaseInMap ( setUp->alpha,
201  setUp->bFactorChange,
202  &maxBandCopy,
203  &maxThetaCopy,
204  &maxPhiCopy,
205  &maxGlIntCopy,
206  setUp,
207  setUp->useCOM,
208  setUp->noIQRsFromMap,
209  setUp->verbose,
210  setUp->clearMapData,
211  setUp->rotChangeDefault,
212  setUp->overlayDefaults,
213  setUp->maskBlurFactor,
214  setUp->maskBlurFactorGiven );
215 
216  inMap->mapPhaselessToSphere ( setUp,
217  maxThetaCopy,
218  maxPhiCopy,
219  minSpacingCopy,
220  setUp->manualShells,
221  true,
222  false );
223 
224  //======================================== Database comparison settings
226 
227  // ... Settings regarding which distances to compute
228  dbCompSetUp->energyLevelDist = true;
229  dbCompSetUp->traceSigmaDist = true;
230  dbCompSetUp->fullRotFnDist = true;
231 
232  // ... Settings regarding hierarchical distance computation
233  dbCompSetUp->enLevelsThreshold = -999.9;
234  dbCompSetUp->trSigmaThreshold = -999.9;
235 
236  // ... Settings regarding the task
237  dbCompSetUp->taskToPerform = ProSHADE::Distances;
238 
239  // ... Settings regarding the database
240  dbCompSetUp->databaseName = setUp->databaseName;
241  dbCompSetUp->volumeTolerance = setUp->volumeTolerance;
242  dbCompSetUp->clearMapData = false;
243  dbCompSetUp->extraSpace = setUp->extraSpace;
244 
245  // ... Settings regarding loudness
246  dbCompSetUp->verbose = setUp->verbose;
247  dbCompSetUp->htmlReport = setUp->htmlReport;
248 
249  //======================================== For each input structure, fragment the map
250  for ( unsigned int strIt = 0; strIt < static_cast<unsigned int> ( strRadii.size ( ) ); strIt++ )
251  {
252  //==================================== Fragmentation set up
253  setUp->mapFragBoxSize = strRadii.at(strIt)[1] * 1.5;
254 
255  //==================================== Do the fragmentation
256  std::vector<ProSHADE_internal::ProSHADE_data*> frags = inMap->fragmentMap ( setUp, false );
257 
258  //==================================== For each passing fragment, check against the input db
259  for ( unsigned int frIt = 0; frIt < static_cast<unsigned int> ( frags.size() ); frIt++ )
260  {
261  //================================ Initialise fragment comparison settings
262  std::stringstream frgName;
263  frgName << "proshade_tmp_frag" << strIt << "_" << frIt << ".map";
264  dbCompSetUp->structFiles.clear ( );
265  dbCompSetUp->structFiles.emplace_back ( frgName.str() );
266 
267  //================================ Print the fragment - this could be made fasted by using the internal representation, but it would also require more work. Will get back to this in due time.
268  frags.at(frIt)->writeMap ( frgName.str(), frags.at(frIt)->getMap ( ) );
269 
270  //================================ Now, compare the fragment against the database
271  ProSHADE::ProSHADE *compareFrag = new ProSHADE::ProSHADE ( dbCompSetUp );
272 
273  //================================ And get distances
274  std::vector<double> crossCorrDists = compareFrag->getCrossCorrDists ( );
275  std::vector<double> traceSigmaDists = compareFrag->getTraceSigmaDists ( );
276  std::vector<double> rotFunDists = compareFrag->getRotFunctionDists ( );
277 
278  //================================ Release the computation object
279  delete compareFrag;
280 
281  for ( int i = 0; i < static_cast<int> ( crossCorrDists.size() ); i++ )
282  {
283 // if ( rotFunDists.at(i) > 0.1 )
284  {
285  printf ( "%i\t%i\t||\t%+.3f\t\t%+.3f\t\t%+.3f\t\t%s\n", strIt, frIt, crossCorrDists.at(i), traceSigmaDists.at(i), rotFunDists.at(i), dbCompSetUp->structFiles.at(i+1).c_str() );
286  }
287  }
288  }
289 
290  std::cout << " ... " << std::endl; exit (0);
291 
292  //==================================== Clean up
293  for ( unsigned int iter = 0; iter < static_cast<unsigned int> ( frags.size ( ) ); iter++ )
294  {
295  delete frags.at(iter);
296  }
297  frags.clear ( );
298  }
299 
300  //======================================== Free database comparison memory
301  delete dbCompSetUp;
302 
303  //======================================== Free memory
304  delete inMap;
305  delete setUp;
306 
307  //======================================== DONE
308  return 0;
309 }
310 
This class deals with reading in the data and computing structure specific information including the ...
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:3019
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
unsigned int checkFileType(std::string fileName)
This function checks the input file for being either PDB or MAP formatted.
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
void mapPhaselessToSphere(ProSHADE::ProSHADE_settings *settings, double theta, double phi, double shellSz, unsigned int manualShells=0, bool keepInMemory=false, bool rotDefaults=false)
This function assumes the data have been processed and maps them onto a set of concentric spheres wit...
bool overlayDefaults
If true, the shell spacing and distances will be doube to their typical values. This is to speed up m...
Definition: ProSHADE.h:176
unsigned int bandwidth
This parameter determines the angular resolution of the spherical harmonics decomposition.
Definition: ProSHADE.h:80
void keepPhaseInMap(double alpha, double bFac, unsigned int *bandwidth, unsigned int *theta, unsigned int *phi, unsigned int *glIntegOrder, ProSHADE::ProSHADE_settings *settings, bool useCom=true, double maxMapIQR=10.0, int verbose=0, bool clearMapData=true, bool rotDefaults=false, bool overlapDefaults=false, double blurFactor=500.0, bool maskBlurFactorGiven=false)
This function keeps the phase information from the density map and prepares the data for SH coefficie...
bool htmlReport
Should HTML report for the run be created?
Definition: ProSHADE.h:185
int verbose
Should the software report on the progress, or just be quiet? Value between 0 (quiet) and 4 (loud) ...
Definition: ProSHADE.h:188
bool saveWithAndWithout
This option decides whether both with and without phase spherical harmonics should be saved...
Definition: ProSHADE.h:102
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
void getDensityMapFromPDB(std::string fileName, double *shellDistance, double resolution, unsigned int *bandwidth, unsigned int *theta, unsigned int *phi, unsigned int *glIntegOrder, double *extraSpace, bool mapResDefault, ProSHADE::ProSHADE_settings *settings, double Bfactor=80.0, bool hpFirstLineCom=false, bool overlayDefaults=false)
Function to read in the PDB file and compute the theoretical density map.
double shellSpacing
This parameter determines how far the radial shells should be from each other.
Definition: ProSHADE.h:96
double volumeTolerance
The percentage tolerance on each dimmension when comparing one structure to entire database...
Definition: ProSHADE.h:159
std::vector< double > getTraceSigmaDists(void)
Accessor function for the trace sigma distances vector.
Definition: ProSHADE.cpp:2993
bool mapResDefault
This variable states if default resolution should be used, or whether the user has supplied a differe...
Definition: ProSHADE.h:85
bool traceSigmaDist
Should the trace sigma distances descriptor be computed.
Definition: ProSHADE.h:133
std::vector< ProSHADE_data * > fragmentMap(ProSHADE::ProSHADE_settings *settings, bool userCOM)
This function takes the map and fragments it into boxes of given size, returning vector of data objec...
bool rotChangeDefault
If map rotation is selected, the default automatic parameter decision is changed. This variable state...
Definition: ProSHADE.h:170
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
void getDensityMapFromMAP(std::string fileName, double *shellDistance, double resolution, unsigned int *bandwidth, unsigned int *theta, unsigned int *phi, unsigned int *glIntegOrder, double *extraSpace, bool mapResDefault, bool rotDefaults, ProSHADE::ProSHADE_settings *settings, bool overlayDefaults=false)
Function to read in the MAP file and provide the basic processing.
double mapFragBoxSize
Should the clear map be fragmented into boxes? If so, put box size here, otherwise leave 0...
Definition: ProSHADE.h:151
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
std::vector< double > getCrossCorrDists(void)
Accessor function for the cross-correlation distances vector.
Definition: ProSHADE.cpp:2967
bool useCOM
Should the Centre of Mass (COM) be used to center the structure in the cell?
Definition: ProSHADE.h:105
double modelRadius(std::string modelPath)
Miscellanous function for obtaining the maximum distance between atoms in a model file...
Definition: ProSHADE.cpp:2931
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
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