ProSHADE  0.7.5.1 (JAN 2021)
Protein Shape Detection
pyProSHADE.cpp
Go to the documentation of this file.
1 
21 //==================================================== Include PyBind11 header
22 #include <pybind11/pybind11.h>
23 #include <pybind11/stl.h>
24 #include <pybind11/numpy.h>
25 
26 //==================================================== Add the ProSHADE_settings and ProSHADE_run classes to the PyBind11 module
27 void add_settingsClass ( pybind11::module& pyProSHADE )
28 {
29  //================================================ Export the ProSHADE_settings class
30  pybind11::class_ < ProSHADE_settings > ( pyProSHADE, "ProSHADE_settings" )
31 
32  //============================================ Constructors (destructors do not need wrappers???)
33  .def ( pybind11::init < > ( ) )
34  .def ( pybind11::init < ProSHADE_Task > ( ), pybind11::arg ( "task" ) )
35 
36  //============================================ Member variables
37  .def_readwrite ( "task", &ProSHADE_settings::task )
38 
39  .def_readwrite ( "inputFiles", &ProSHADE_settings::inputFiles )
40  .def_readwrite ( "forceP1", &ProSHADE_settings::forceP1 )
41  .def_readwrite ( "removeWaters", &ProSHADE_settings::removeWaters )
42  .def_readwrite ( "firstModelOnly", &ProSHADE_settings::firstModelOnly )
43 
44  .def_readwrite ( "requestedResolution", &ProSHADE_settings::requestedResolution )
45  .def_readwrite ( "changeMapResolution", &ProSHADE_settings::changeMapResolution )
46  .def_readwrite ( "changeMapResolutionTriLinear", &ProSHADE_settings::changeMapResolutionTriLinear )
47 
48  .def_readwrite ( "pdbBFactorNewVal", &ProSHADE_settings::pdbBFactorNewVal )
49 
50  .def_readwrite ( "maxBandwidth", &ProSHADE_settings::maxBandwidth )
51  .def_readwrite ( "rotationUncertainty", &ProSHADE_settings::rotationUncertainty )
52 
53  .def_readwrite ( "usePhase", &ProSHADE_settings::usePhase )
54 
55  .def_readwrite ( "maxSphereDists", &ProSHADE_settings::maxSphereDists )
56 
57  .def_readwrite ( "integOrder", &ProSHADE_settings::integOrder )
58  .def_readwrite ( "taylorSeriesCap", &ProSHADE_settings::taylorSeriesCap )
59 
60  .def_readwrite ( "normaliseMap", &ProSHADE_settings::normaliseMap )
61 
62  .def_readwrite ( "invertMap", &ProSHADE_settings::invertMap )
63 
64  .def_readwrite ( "blurFactor", &ProSHADE_settings::blurFactor )
65  .def_readwrite ( "maskingThresholdIQRs", &ProSHADE_settings::maskingThresholdIQRs )
66  .def_readwrite ( "maskMap", &ProSHADE_settings::maskMap )
67  .def_readwrite ( "useCorrelationMasking", &ProSHADE_settings::useCorrelationMasking )
68  .def_readwrite ( "halfMapKernel", &ProSHADE_settings::halfMapKernel )
69  .def_readwrite ( "correlationKernel", &ProSHADE_settings::correlationKernel )
70  .def_readwrite ( "saveMask", &ProSHADE_settings::saveMask )
71  .def_readwrite ( "maskFileName", &ProSHADE_settings::maskFileName )
72 
73  .def_readwrite ( "reBoxMap", &ProSHADE_settings::reBoxMap )
74  .def_readwrite ( "boundsExtraSpace", &ProSHADE_settings::boundsExtraSpace )
75  .def_readwrite ( "boundsSimilarityThreshold", &ProSHADE_settings::boundsSimilarityThreshold )
76  .def_readwrite ( "useSameBounds", &ProSHADE_settings::useSameBounds )
77  .def_readwrite ( "forceBounds", &ProSHADE_settings::forceBounds )
78 
79  .def_readwrite ( "moveToCOM", &ProSHADE_settings::moveToCOM )
80 
81  .def_readwrite ( "addExtraSpace", &ProSHADE_settings::addExtraSpace )
82 
83  .def_readwrite ( "progressiveSphereMapping", &ProSHADE_settings::progressiveSphereMapping )
84 
85  .def_readwrite ( "outName", &ProSHADE_settings::outName )
86 
87  .def_readwrite ( "computeEnergyLevelsDesc", &ProSHADE_settings::computeEnergyLevelsDesc )
88  .def_readwrite ( "enLevMatrixPowerWeight", &ProSHADE_settings::enLevMatrixPowerWeight )
89  .def_readwrite ( "computeTraceSigmaDesc", &ProSHADE_settings::computeTraceSigmaDesc )
90  .def_readwrite ( "computeRotationFuncDesc", &ProSHADE_settings::computeRotationFuncDesc )
91 
92  .def_readwrite ( "peakNeighbours", &ProSHADE_settings::peakNeighbours )
93  .def_readwrite ( "noIQRsFromMedianNaivePeak", &ProSHADE_settings::noIQRsFromMedianNaivePeak )
94 
95  .def_readwrite ( "smoothingFactor", &ProSHADE_settings::smoothingFactor )
96 
97  .def_readwrite ( "symMissPeakThres", &ProSHADE_settings::symMissPeakThres )
98  .def_readwrite ( "axisErrTolerance", &ProSHADE_settings::axisErrTolerance )
99  .def_readwrite ( "axisErrToleranceDefault", &ProSHADE_settings::axisErrToleranceDefault )
100  .def_readwrite ( "minSymPeak", &ProSHADE_settings::minSymPeak )
101  .def_readwrite ( "recommendedSymmetryType", &ProSHADE_settings::recommendedSymmetryType )
102  .def_readwrite ( "recommendedSymmetryFold", &ProSHADE_settings::recommendedSymmetryFold )
103  .def_readwrite ( "requestedSymmetryType", &ProSHADE_settings::requestedSymmetryType )
104  .def_readwrite ( "requestedSymmetryFold", &ProSHADE_settings::requestedSymmetryFold )
105  .def_readwrite ( "usePeakSearchInRotationFunctionSpace", &ProSHADE_settings::usePeakSearchInRotationFunctionSpace)
106  .def_readwrite ( "useBiCubicInterpolationOnPeaks", &ProSHADE_settings::useBiCubicInterpolationOnPeaks )
107  .def_readwrite ( "maxSymmetryFold", &ProSHADE_settings::maxSymmetryFold )
108 
109  .def_readwrite ( "overlayStructureName", &ProSHADE_settings::overlayStructureName )
110  .def_readwrite ( "rotTrsJSONFile", &ProSHADE_settings::rotTrsJSONFile )
111 
112  .def_readwrite ( "verbose", &ProSHADE_settings::verbose )
113 
114  //============================================ Mutators
115  .def ( "addStructure", &ProSHADE_settings::addStructure, "Adds a structure file name to the appropriate variable.", pybind11::arg ( "structure" ) )
116  .def ( "setResolution", &ProSHADE_settings::setResolution, "This function sets the resolution in the appropriate variable.", pybind11::arg ( "resolution" ) )
117  .def ( "setPDBBFactor", &ProSHADE_settings::setPDBBFactor, "Sets the requested B-factor value for PDB files in the appropriate variable.", pybind11::arg ( "newBF" ) )
118  .def ( "setNormalisation", &ProSHADE_settings::setNormalisation, "Sets the requested map normalisation value in the appropriate variable.", pybind11::arg ( "normalise" ) )
119  .def ( "setMapInversion", &ProSHADE_settings::setMapInversion, "Sets the requested map inversion value in the appropriate variable.", pybind11::arg ( "mInv" ) )
120  .def ( "setVerbosity", &ProSHADE_settings::setVerbosity, "Sets the requested verbosity in the appropriate variable.", pybind11::arg ( "verbosity" ) )
121  .def ( "setMaskBlurFactor", &ProSHADE_settings::setMaskBlurFactor, "Sets the requested map blurring factor in the appropriate variable.", pybind11::arg ( "blurFac" ) )
122  .def ( "setMaskIQR", &ProSHADE_settings::setMaskIQR, "Sets the requested number of IQRs for masking threshold in the appropriate variable.", pybind11::arg ( "noIQRs" ) )
123  .def ( "setMasking", &ProSHADE_settings::setMasking, "Sets the requested map masking decision value in the appropriate variable.", pybind11::arg ( "mask" ) )
124  .def ( "setCorrelationMasking", &ProSHADE_settings::setCorrelationMasking, "Sets the requested map masking type in the appropriate variable.", pybind11::arg ( "corMask" ) )
125  .def ( "setTypicalNoiseSize", &ProSHADE_settings::setTypicalNoiseSize, "Sets the requested \"fake\" half-map kernel in the appropriate variable.", pybind11::arg ( "typNoi" ) )
126  .def ( "setMinimumMaskSize", &ProSHADE_settings::setMinimumMaskSize, "Sets the requested minimum mask size.", pybind11::arg ( "minMS" ) )
127  .def ( "setMaskSaving", &ProSHADE_settings::setMaskSaving, "Sets whether the mask should be saved.", pybind11::arg ( "savMsk" ) )
128  .def ( "setMaskFilename", &ProSHADE_settings::setMaskFilename, "Sets where the mask should be saved.", pybind11::arg ( "mskFln" ) )
129  .def ( "setMapReboxing", &ProSHADE_settings::setMapReboxing, "Sets whether re-boxing needs to be done in the appropriate variable.", pybind11::arg ( "reBx" ) )
130  .def ( "setBoundsSpace", &ProSHADE_settings::setBoundsSpace, "Sets the requested number of angstroms for extra space in re-boxing in the appropriate variable.", pybind11::arg ( "boundsExSp" ) )
131  .def ( "setBoundsThreshold", &ProSHADE_settings::setBoundsThreshold, "Sets the threshold for number of indices difference acceptable to make index sizes same in the appropriate variable.", pybind11::arg ( "boundsThres" ) )
132  .def ( "setSameBoundaries", &ProSHADE_settings::setSameBoundaries, "Sets whether same boundaries should be used in the appropriate variable.", pybind11::arg ( "sameB" ) )
133  .def ( "setOutputFilename", &ProSHADE_settings::setOutputFilename, "Sets the requested output file name in the appropriate variable.", pybind11::arg ( "oFileName" ) )
134  .def ( "setMapResolutionChange", &ProSHADE_settings::setMapResolutionChange, "Sets the requested map resolution change decision in the appropriate variable.", pybind11::arg ( "mrChange" ) )
135  .def ( "setMapResolutionChangeTriLinear", &ProSHADE_settings::setMapResolutionChangeTriLinear, "Sets the requested map resolution change decision using tri-linear interpolation in the appropriate variable.", pybind11::arg ( "mrChange" ) )
136  .def ( "setMapCentering", &ProSHADE_settings::setMapCentering, "Sets the requested map centering decision value in the appropriate variable.", pybind11::arg ( "com" ) )
137  .def ( "setExtraSpace", &ProSHADE_settings::setExtraSpace, "Sets the requested map extra space value in the appropriate variable.", pybind11::arg ( "exSpace" ) )
138  .def ( "setBandwidth", &ProSHADE_settings::setBandwidth, "Sets the requested spherical harmonics bandwidth in the appropriate variable.", pybind11::arg ( "band" ) )
139  .def ( "setProgressiveSphereMapping", &ProSHADE_settings::setProgressiveSphereMapping, "Sets the requested sphere mapping value settings approach in the appropriate variable.", pybind11::arg ( "progSphMap" ) )
140  .def ( "setSphereDistances", &ProSHADE_settings::setSphereDistances, "Sets the requested distance between spheres in the appropriate variable.", pybind11::arg ( "sphDist" ) )
141  .def ( "setIntegrationOrder", &ProSHADE_settings::setIntegrationOrder, "Sets the requested order for the Gauss-Legendre integration in the appropriate variable.", pybind11::arg ( "intOrd" ) )
142  .def ( "setTaylorSeriesCap", &ProSHADE_settings::setTaylorSeriesCap, "Sets the requested Taylor series cap for the Gauss-Legendre integration in the appropriate variable.", pybind11::arg ( "tayCap" ) )
143  .def ( "setEnergyLevelsComputation", &ProSHADE_settings::setEnergyLevelsComputation, "Sets whether the energy level distance descriptor should be computed.", pybind11::arg ( "enLevDesc" ) )
144  .def ( "setTraceSigmaComputation", &ProSHADE_settings::setTraceSigmaComputation, "Sets whether the trace sigma distance descriptor should be computed.", pybind11::arg ( "trSigVal" ) )
145  .def ( "setRotationFunctionComputation", &ProSHADE_settings::setRotationFunctionComputation, "Sets whether the rotation function distance descriptor should be computed.", pybind11::arg ( "rotfVal" ) )
146  .def ( "setPeakNeighboursNumber", &ProSHADE_settings::setPeakNeighboursNumber, "Sets the number of neighbour values that have to be smaller for an index to be considered a peak.", pybind11::arg ( "pkS" ) )
147  .def ( "setPeakNaiveNoIQR", &ProSHADE_settings::setPeakNaiveNoIQR, "Sets the number of IQRs from the median for threshold height a peak needs to be considered a peak.", pybind11::arg ( "noIQRs" ) )
148  .def ( "setPhaseUsage", &ProSHADE_settings::setPhaseUsage, "Sets whether the phase information will be used.", pybind11::arg ( "phaseUsage" ) )
149  .def ( "setEnLevShellWeight", &ProSHADE_settings::setEnLevShellWeight, "Sets the weight of shell position for the energy levels computation.", pybind11::arg ( "mPower" ) )
150  .def ( "setGroupingSmoothingFactor", &ProSHADE_settings::setGroupingSmoothingFactor, "Sets the grouping smoothing factor into the proper variable.", pybind11::arg ( "smFact" ) )
151  .def ( "setMissingPeakThreshold", &ProSHADE_settings::setMissingPeakThreshold, "Sets the threshold for starting the missing peaks procedure.", pybind11::arg ( "mpThres" ) )
152  .def ( "setAxisComparisonThreshold", &ProSHADE_settings::setAxisComparisonThreshold, "Sets the threshold for matching symmetry axes.", pybind11::arg ( "axThres" ) )
153  .def ( "setAxisComparisonThresholdBehaviour", &ProSHADE_settings::setAxisComparisonThresholdBehaviour, "Sets the automatic symmetry axis tolerance decreasing.", pybind11::arg ( "behav" ) )
154  .def ( "setMinimumPeakForAxis", &ProSHADE_settings::setMinimumPeakForAxis, "Sets the minimum peak height for symmetry axis to be considered.", pybind11::arg ( "minSP" ) )
155  .def ( "setRecommendedSymmetry", &ProSHADE_settings::setRecommendedSymmetry, "Sets the ProSHADE detected symmetry type.", pybind11::arg ( "val" ) )
156  .def ( "setRecommendedFold", &ProSHADE_settings::setRecommendedFold, "Sets the ProSHADE detected symmetry fold.", pybind11::arg ( "val" ) )
157  .def ( "setRequestedSymmetry", &ProSHADE_settings::setRequestedSymmetry, "Sets the user requested symmetry type.", pybind11::arg ( "val" ) )
158  .def ( "setRequestedFold", &ProSHADE_settings::setRequestedFold, "Sets the user requested symmetry fold.", pybind11::arg ( "val" ) )
159  .def ( "setDetectedSymmetry", &ProSHADE_settings::setDetectedSymmetry, "Sets the final detected symmetry axes information.", pybind11::arg ( "sym" ) )
160  .def ( "setOverlaySaveFile", &ProSHADE_settings::setOverlaySaveFile, "Sets the filename to which the overlay structure is to be save into.", pybind11::arg ( "filename" ) )
161  .def ( "setOverlayJsonFile", &ProSHADE_settings::setOverlayJsonFile, "Sets the filename to which the overlay operations are to be save into.", pybind11::arg ( "filename" ) )
162  .def ( "setSymmetryRotFunPeaks", &ProSHADE_settings::setSymmetryRotFunPeaks, "Sets the symmetry detection algorithm type.", pybind11::arg ( "rotFunPeaks" ) )
163  .def ( "setBicubicInterpolationSearch", &ProSHADE_settings::setBicubicInterpolationSearch, "Sets the bicubic interpolation on peaks.", pybind11::arg ( "bicubPeaks" ) )
164  .def ( "setMaxSymmetryFold", &ProSHADE_settings::setMaxSymmetryFold, "Sets the maximum symmetry fold (well, the maximum prime symmetry fold).", pybind11::arg ( "maxFold" ) )
165 
166  //============================================ Command line parsing
167  .def ( "getCommandLineParams",
168  [] ( ProSHADE_settings &self, std::vector < std::string > args )
169  {
170  std::vector < char * > cstrs; cstrs.reserve ( args.size() );
171 
172  for ( auto &s : args )
173  cstrs.push_back ( const_cast < char * > ( s.c_str ( ) ) );
174 
175  return self.getCommandLineParams ( cstrs.size ( ), cstrs.data ( ) );
176  }, "This function takes a VectorOfStrings and parses it as if it were command line arguments, filling in the calling ProSHADE_settings class with the values." )
177 
178  //============================================ Debugging
179  .def ( "printSettings", &ProSHADE_settings::printSettings, "This function prints the current values in the settings object." )
180 
181  //============================================ Description
182  .def ( "__repr__", [] ( const ProSHADE_settings &a ) { return "<ProSHADE_settings class object> (Settings class is used to set all settings values in a single place)"; } );
183 
184  //================================================ Export the ProSHADE_run class
185  pybind11::class_ < ProSHADE_run > ( pyProSHADE, "ProSHADE_run" )
186 
187  //============================================ Constructors (destructors do not need wrappers???)
188  .def ( pybind11::init < ProSHADE_settings* > ( ) )
189 
190  //============================================ General accessors
191  .def ( "getNoStructures", &ProSHADE_run::getNoStructures, "This function returns the number of structures used." )
192  .def ( "getVerbose", &ProSHADE_run::getVerbose, "This function returns the verbose value." )
193 
194  //============================================ Distances results accessor functions wrapped as lambda functions for numpy return types
195  .def ( "getEnergyLevelsVector",
196  [] ( ProSHADE_run &self ) -> pybind11::array_t < float >
197  {
198  //== Get the values
199  std::vector< proshade_double > vals = self.getEnergyLevelsVector ();
200 
201  //== Allocate memory for the numpy values
202  float* npVals = new float[static_cast<unsigned int> (vals.size())];
203  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
204 
205  //== Copy values
206  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = vals.at(iter); }
207 
208  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
209  pybind11::capsule pyCapsuleEnLevs ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
210 
211  //== Copy the value
212  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<unsigned int> (vals.size()) }, // Shape
213  { sizeof(float) }, // C-stype strides
214  npVals, // Data
215  pyCapsuleEnLevs ); // Capsule (C++ destructor, basically)
216 
217  //== Done
218  return ( retArr );
219  }, "This function returns the energy level distances vector from the first to all other structures." )
220 
221  .def ( "getTraceSigmaVector",
222  [] ( ProSHADE_run &self ) -> pybind11::array_t < float >
223  {
224  //== Get the values
225  std::vector< proshade_double > vals = self.getTraceSigmaVector ();
226 
227  //== Allocate memory for the numpy values
228  float* npVals = new float[static_cast<unsigned int> (vals.size())];
229  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
230 
231  //== Copy values
232  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = vals.at(iter); }
233 
234  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
235  pybind11::capsule pyCapsuleTrSigs ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
236 
237  //== Copy the value
238  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<unsigned int> (vals.size()) }, // Shape
239  { sizeof(float) }, // C-stype strides
240  npVals, // Data
241  pyCapsuleTrSigs ); // Capsule (C++ destructor, basically)
242 
243  //== Done
244  return ( retArr );
245  }, "This function returns the trace sigma distances vector from the first to all other structures." )
246 
247  .def ( "getRotationFunctionVector",
248  [] ( ProSHADE_run &self ) -> pybind11::array_t < float >
249  {
250  //== Get the values
251  std::vector< proshade_double > vals = self.getRotationFunctionVector ();
252 
253  //== Allocate memory for the numpy values
254  float* npVals = new float[static_cast<unsigned int> (vals.size())];
255  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
256 
257  //== Copy values
258  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = vals.at(iter); }
259 
260  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
261  pybind11::capsule pyCapsuleRotFun ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
262 
263  //== Copy the value
264  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<unsigned int> (vals.size()) }, // Shape
265  { sizeof(float) }, // C-stype strides
266  npVals, // Data
267  pyCapsuleRotFun ); // Capsule (C++ destructor, basically)
268 
269  //== Done
270  return ( retArr );
271  }, "This function returns the full rotation function distances vector from the first to all other structures." )
272 
273  //============================================ Symmetry results accessor functions
274  .def ( "getSymmetryType", &ProSHADE_run::getSymmetryType, "This is the main accessor function for the user to get to know what symmetry type ProSHADE has detected and recommends." )
275  .def ( "getSymmetryFold", &ProSHADE_run::getSymmetryFold, "This is the main accessor function for the user to get to know what symmetry fold ProSHADE has detected and recommends." )
276  .def ( "getSymmetryAxis", &ProSHADE_run::getSymmetryAxis, "This function returns a single symmetry axis as a vector of strings from the recommended symmetry axes list.", pybind11::arg ( "axisNo" ) )
277  .def ( "getAllCSyms",
278  [] ( ProSHADE_run &self ) -> pybind11::array_t < float >
279  {
280  //== Get the values
281  std::vector< std::vector< proshade_double > > vals = self.getAllCSyms ();
282 
283  //== Allocate memory for the numpy values
284  float* npVals = new float[static_cast<unsigned int> ( vals.size() * 6 )];
285  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
286 
287  //== Copy values
288  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { for ( proshade_unsign it = 0; it < 6; it++ ) { npVals[(iter*6)+it] = vals.at(iter).at(it); } }
289 
290  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
291  pybind11::capsule pyCapsuleSymList ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
292 
293  //== Copy the value
294  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<int> ( vals.size() ), static_cast<int> ( 6 ) }, // Shape
295  { 6 * sizeof(float), sizeof(float) }, // C-stype strides
296  npVals, // Data
297  pyCapsuleSymList ); // Capsule
298 
299  //== Done
300  return ( retArr );
301  }, "This function returns a all symmetry axes as a 2D numpy array." )
302 
303  //============================================ Reboxing results accessor functions as lambda functions directly returning numpy arrays
304  .def ( "getOriginalBounds",
305  [] ( ProSHADE_run &self, proshade_unsign strNo ) -> pybind11::array_t < float >
306  {
307  //== Get values
308  std::vector< proshade_signed > vals = self.getOriginalBounds ( strNo );
309 
310  //== Allocate memory for the numpy values
311  float* npVals = new float[static_cast<proshade_unsign> ( vals.size() )];
312  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
313 
314  //== Copy values
315  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = vals.at(iter); }
316 
317  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
318  pybind11::capsule pyCapsuleOrigBnds ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
319 
320  //== Copy the value
321  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<proshade_unsign> ( vals.size() ) }, // Shape
322  { sizeof(float) }, // C-stype strides
323  npVals, // Data
324  pyCapsuleOrigBnds ); // Capsule (C++ destructor, basically)
325 
326  //== Done
327  return ( retArr );
328  }, "This function returns the original structure boundaries as numpy array." )
329 
330  .def ( "getReBoxedBounds",
331  [] ( ProSHADE_run &self, proshade_unsign strNo ) -> pybind11::array_t < float >
332  {
333  //== Get values
334  std::vector< proshade_signed > vals = self.getReBoxedBounds ( strNo );
335 
336  //== Allocate memory for the numpy values
337  float* npVals = new float[static_cast<proshade_unsign> ( vals.size() )];
338  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
339 
340  //== Copy values
341  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = vals.at(iter); }
342 
343  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
344  pybind11::capsule pyCapsuleReBoBnds ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
345 
346  //== Copy the value
347  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<proshade_unsign> ( vals.size() ) }, // Shape
348  { sizeof(float) }, // C-stype strides
349  npVals, // Data
350  pyCapsuleReBoBnds ); // Capsule (C++ destructor, basically)
351 
352  //== Done
353  return ( retArr );
354  }, "This function returns the re-boxed structure boundaries as numpy array." )
355 
356 
357  .def ( "getReBoxedMap",
358  [] ( ProSHADE_run &self, proshade_unsign strNo ) -> pybind11::array_t < float >
359  {
360  //== Get the values
361  std::vector< proshade_signed > vals = self.getReBoxedBounds ( strNo );
362 
363  //== Determine dimensions
364  proshade_unsign xDim = vals.at(1) - vals.at(0) + 1;
365  proshade_unsign yDim = vals.at(3) - vals.at(2) + 1;
366  proshade_unsign zDim = vals.at(5) - vals.at(4) + 1;
367 
368  //== Allocate memory for the numpy values
369  float* npVals = new float[xDim * yDim * zDim];
370  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
371 
372  //== Copy values
373  for ( proshade_unsign iter = 0; iter < (xDim * yDim * zDim); iter++ ) { npVals[iter] = self.getMapValue ( strNo, iter ); }
374 
375  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
376  pybind11::capsule pyCapsuleRebMap ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
377 
378  //== Copy the value
379  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { xDim, yDim, zDim }, // Shape
380  { yDim * zDim * sizeof(float), zDim * sizeof(float), sizeof(float) }, // C-stype strides
381  npVals, // Data
382  pyCapsuleRebMap ); // Capsule
383 
384  //== Done
385  return ( retArr );
386  }, "This function returns the re-boxed structure map as a numpy 3D array." )
387 
388  //============================================ Overlay results accessor functions as lambda functions directly returning numpy arrays
389  .def ( "getEulerAngles",
390  [] ( ProSHADE_run &self ) -> pybind11::array_t < float >
391  {
392  //== Get the values
393  std::vector< proshade_double > vals = self.getEulerAngles ( );
394 
395  //== Allocate memory for the numpy values
396  float* npVals = new float[static_cast<proshade_unsign> ( vals.size() )];
397  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
398 
399  //== Copy values
400  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = vals.at(iter); }
401 
402  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
403  pybind11::capsule pyCapsuleEulAngs ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
404 
405  //== Copy the value
406  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<proshade_unsign> ( vals.size() ) }, // Shape
407  { sizeof(float) }, // C-stype strides
408  npVals, // Data
409  pyCapsuleEulAngs ); // Capsule
410 
411  //== Done
412  return ( retArr );
413  }, "This function returns the vector of Euler angles with best overlay correlation." )
414 
415  .def ( "getOptimalRotMat",
416  [] ( ProSHADE_run &self ) -> pybind11::array_t < float >
417  {
418  //== Get the values
419  std::vector< proshade_double > vals = self.getOptimalRotMat ( );
420 
421  //== Allocate memory for the numpy values
422  float* npVals = new float[static_cast<proshade_unsign> ( vals.size() )];
423  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
424 
425  //== Copy values
426  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = vals.at(iter); }
427 
428  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
429  pybind11::capsule pyCapsuleRotMat ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
430 
431  //== Copy the value
432  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { 3, 3 }, // Shape
433  { 3 * sizeof(float), sizeof(float) }, // C-stype strides
434  npVals, // Data
435  pyCapsuleRotMat ); // Capsule
436 
437  //== Done
438  return ( retArr );
439  }, "This function returns the vector of Euler angles with best overlay correlation." )
440 
441  .def ( "getTranslationToOrigin",
442  [] ( ProSHADE_run &self ) -> pybind11::array_t < float >
443  {
444  //== Get the values
445  std::vector< proshade_double > vals = self.getTranslationToOrigin ( );
446 
447  //== Allocate memory for the numpy values
448  float* npVals = new float[static_cast<proshade_unsign> ( vals.size() )];
449  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
450 
451  //== Copy values
452  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = vals.at(iter); }
453 
454  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
455  pybind11::capsule pyCapsuleTTO ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
456 
457  //== Copy the value
458  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<proshade_unsign> ( vals.size() ) }, // Shape
459  { sizeof(float) }, // C-stype strides
460  npVals, // Data
461  pyCapsuleTTO ); // Capsule
462 
463  //== Done
464  return ( retArr );
465  }, "This function returns the negative values of the position of the rotation centre (the point about which the rotation should be done)." )
466  .def ( "getOriginToOverlayTranslation",
467  [] ( ProSHADE_run &self ) -> pybind11::array_t < float >
468  {
469  //== Get the values
470  std::vector< proshade_double > vals = self.getOriginToOverlayTranslation ( );
471 
472  //== Allocate memory for the numpy values
473  float* npVals = new float[static_cast<proshade_unsign> ( vals.size() )];
474  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
475 
476  //== Copy values
477  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = vals.at(iter); }
478 
479  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
480  pybind11::capsule pyCapsuleOTOT ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
481 
482  //== Copy the value
483  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<proshade_unsign> ( vals.size() ) }, // Shape
484  { sizeof(float) }, // C-stype strides
485  npVals, // Data
486  pyCapsuleOTOT ); // Capsule
487 
488  //== Done
489  return ( retArr );
490  }, "This function returns the translation required to move the structure from origin to optimal overlay." )
491 
492 
493  //============================================ Description
494  .def ( "__repr__", [] ( const ProSHADE_run &a ) { return "<ProSHADE_run class object> (Run class constructor takes a ProSHADE_settings object and completes a single run according to the settings object information)"; } );
495 }
ProSHADE_settings::noIQRsFromMedianNaivePeak
proshade_double noIQRsFromMedianNaivePeak
When doing peak searching, how many IQRs from the median the threshold for peak height should be (in ...
Definition: ProSHADE_settings.hpp:162
ProSHADE_settings::setOverlayJsonFile
void setOverlayJsonFile(std::string filename)
Sets the filename to which the overlay operations are to be save into.
Definition: ProSHADE.cpp:1120
ProSHADE_settings::maxBandwidth
proshade_unsign maxBandwidth
The bandwidth of spherical harmonics decomposition for the largest sphere.
Definition: ProSHADE_settings.hpp:106
ProSHADE_settings::integOrder
proshade_unsign integOrder
The order required for full Gauss-Legendre integration between the spheres.
Definition: ProSHADE_settings.hpp:116
ProSHADE_settings::recommendedSymmetryType
std::string recommendedSymmetryType
The symmetry type that ProSHADE finds the best fitting for the structure. Possible values are "" for ...
Definition: ProSHADE_settings.hpp:172
ProSHADE_settings::rotTrsJSONFile
std::string rotTrsJSONFile
The filename to which the rotation and translation operations are to be saved into.
Definition: ProSHADE_settings.hpp:182
ProSHADE_settings::setEnLevShellWeight
void setEnLevShellWeight(proshade_double mPower)
Sets the weight of shell position for the energy levels computation.
Definition: ProSHADE.cpp:903
ProSHADE_settings::computeTraceSigmaDesc
bool computeTraceSigmaDesc
If true, the trace sigma descriptor will be computed, otherwise all its computations will be omitted.
Definition: ProSHADE_settings.hpp:157
ProSHADE_settings::setTraceSigmaComputation
void setTraceSigmaComputation(bool trSigVal)
Sets whether the trace sigma distance descriptor should be computed.
Definition: ProSHADE.cpp:815
ProSHADE_settings::computeRotationFuncDesc
bool computeRotationFuncDesc
If true, the rotation function descriptor will be computed, otherwise all its computations will be om...
Definition: ProSHADE_settings.hpp:158
ProSHADE_settings::setSameBoundaries
void setSameBoundaries(bool sameB)
Sets whether same boundaries should be used in the appropriate variable.
Definition: ProSHADE.cpp:619
ProSHADE_settings::maxSymmetryFold
proshade_unsign maxSymmetryFold
The highest symmetry fold to search for.
Definition: ProSHADE_settings.hpp:178
ProSHADE_settings::taylorSeriesCap
proshade_unsign taylorSeriesCap
The max limit on the Taylor series expansion done for the abscissas of the Gauss-Legendre integration...
Definition: ProSHADE_settings.hpp:117
ProSHADE_settings::setMinimumMaskSize
void setMinimumMaskSize(proshade_single minMS)
Sets the requested minimum mask size.
Definition: ProSHADE.cpp:521
ProSHADE_settings::boundsExtraSpace
proshade_single boundsExtraSpace
The number of extra angstroms to be added to all re-boxing bounds just for safety.
Definition: ProSHADE_settings.hpp:137
ProSHADE_settings::setBicubicInterpolationSearch
void setBicubicInterpolationSearch(bool bicubPeaks)
Sets the bicubic interpolation on peaks.
Definition: ProSHADE.cpp:1148
ProSHADE_settings::setPeakNaiveNoIQR
void setPeakNaiveNoIQR(proshade_double noIQRs)
Sets the number of IQRs from the median for threshold height a peak needs to be considered a peak.
Definition: ProSHADE.cpp:867
ProSHADE_settings::outName
std::string outName
The file name where the output structure(s) should be saved.
Definition: ProSHADE_settings.hpp:152
ProSHADE_settings::setMapInversion
void setMapInversion(bool mInv)
Sets the requested map inversion value in the appropriate variable.
Definition: ProSHADE.cpp:405
ProSHADE_settings::setMapResolutionChange
void setMapResolutionChange(bool mrChange)
Sets the requested map resolution change decision in the appropriate variable.
Definition: ProSHADE.cpp:652
ProSHADE_settings::setPeakNeighboursNumber
void setPeakNeighboursNumber(proshade_unsign pkS)
Sets the number of neighbour values that have to be smaller for an index to be considered a peak.
Definition: ProSHADE.cpp:849
ProSHADE_settings::blurFactor
proshade_single blurFactor
This is the amount by which B-factors should be increased to create the blurred map for masking.
Definition: ProSHADE_settings.hpp:126
ProSHADE_settings::maskMap
bool maskMap
Should the map be masked from noise?
Definition: ProSHADE_settings.hpp:128
ProSHADE_settings::requestedSymmetryFold
proshade_unsign requestedSymmetryFold
The fold of the requested symmetry (only applicable to C and D symmetry types).
Definition: ProSHADE_settings.hpp:175
ProSHADE_settings::requestedSymmetryType
std::string requestedSymmetryType
The symmetry type requested by the user. Allowed values are C, D, T, O and I.
Definition: ProSHADE_settings.hpp:174
ProSHADE_settings::setSphereDistances
void setSphereDistances(proshade_single sphDist)
Sets the requested distance between spheres in the appropriate variable.
Definition: ProSHADE.cpp:748
ProSHADE_settings::setVerbosity
void setVerbosity(proshade_signed verbosity)
Sets the requested verbosity in the appropriate variable.
Definition: ProSHADE.cpp:421
ProSHADE_settings::setMinimumPeakForAxis
void setMinimumPeakForAxis(proshade_double minSP)
Sets the minimum peak height for symmetry axis to be considered.
Definition: ProSHADE.cpp:991
ProSHADE_settings::saveMask
bool saveMask
Should the mask be saved?
Definition: ProSHADE_settings.hpp:132
ProSHADE_settings::requestedResolution
proshade_single requestedResolution
The resolution to which the calculations are to be done.
Definition: ProSHADE_settings.hpp:98
ProSHADE_settings::minSymPeak
proshade_double minSymPeak
Minimum average peak for symmetry axis to be considered as "real".
Definition: ProSHADE_settings.hpp:171
ProSHADE_settings::setAxisComparisonThresholdBehaviour
void setAxisComparisonThresholdBehaviour(bool behav)
Sets the automatic symmetry axis tolerance decreasing.
Definition: ProSHADE.cpp:974
ProSHADE_settings::setDetectedSymmetry
void setDetectedSymmetry(proshade_double *sym)
Sets the final detected symmetry axes information.
Definition: ProSHADE.cpp:1077
ProSHADE_settings::progressiveSphereMapping
bool progressiveSphereMapping
If true, each shell will have its own angular resolution dependent on the actual number of map points...
Definition: ProSHADE_settings.hpp:149
ProSHADE_settings::maxSphereDists
proshade_single maxSphereDists
The distance between spheres in spherical mapping for the largest sphere.
Definition: ProSHADE_settings.hpp:113
ProSHADE_settings::symMissPeakThres
proshade_double symMissPeakThres
Percentage of peaks that could be missing that would warrant starting the missing peaks search proced...
Definition: ProSHADE_settings.hpp:168
ProSHADE_settings::setPDBBFactor
void setPDBBFactor(proshade_double newBF)
Sets the requested B-factor value for PDB files in the appropriate variable.
Definition: ProSHADE.cpp:373
ProSHADE_settings::addStructure
void addStructure(std::string structure)
Adds a structure file name to the appropriate variable.
Definition: ProSHADE.cpp:341
ProSHADE_settings::setMaskIQR
void setMaskIQR(proshade_single noIQRs)
Sets the requested number of IQRs for masking threshold in the appropriate variable.
Definition: ProSHADE.cpp:454
ProSHADE_settings::setMissingPeakThreshold
void setMissingPeakThreshold(proshade_double mpThres)
Sets the threshold for starting the missing peaks procedure.
Definition: ProSHADE.cpp:938
ProSHADE_run::getVerbose
proshade_signed getVerbose(void)
This function returns the verbose value.
Definition: ProSHADE.cpp:2330
ProSHADE_settings::maskingThresholdIQRs
proshade_single maskingThresholdIQRs
Number of inter-quartile ranges from the median to be used for thresholding the blurred map for maski...
Definition: ProSHADE_settings.hpp:127
ProSHADE_settings::useCorrelationMasking
bool useCorrelationMasking
Should the blurring masking (false) or the correlation masking (true) be used?
Definition: ProSHADE_settings.hpp:129
ProSHADE_settings::usePhase
bool usePhase
If true, the full data will be used, if false, Patterson maps will be used instead and phased data wi...
Definition: ProSHADE_settings.hpp:110
ProSHADE_settings::setMapCentering
void setMapCentering(bool com)
Sets the requested map centering decision value in the appropriate variable.
Definition: ProSHADE.cpp:684
ProSHADE_settings::setMaskBlurFactor
void setMaskBlurFactor(proshade_single blurFac)
Sets the requested map blurring factor in the appropriate variable.
Definition: ProSHADE.cpp:437
ProSHADE_settings::setPhaseUsage
void setPhaseUsage(bool phaseUsage)
Sets whether the phase information will be used.
Definition: ProSHADE.cpp:885
ProSHADE_settings::maskFileName
std::string maskFileName
The filename to which mask should be saved.
Definition: ProSHADE_settings.hpp:133
ProSHADE_settings::verbose
proshade_signed verbose
Should the software report on the progress, or just be quiet? Value between -1 (nothing) and 4 (loud)
Definition: ProSHADE_settings.hpp:185
ProSHADE_settings::setBandwidth
void setBandwidth(proshade_unsign band)
Sets the requested spherical harmonics bandwidth in the appropriate variable.
Definition: ProSHADE.cpp:732
ProSHADE_settings::setBoundsSpace
void setBoundsSpace(proshade_single boundsExSp)
Sets the requested number of angstroms for extra space in re-boxing in the appropriate variable.
Definition: ProSHADE.cpp:586
ProSHADE_settings::changeMapResolutionTriLinear
bool changeMapResolutionTriLinear
Should maps be re-sampled to obtain the required resolution?
Definition: ProSHADE_settings.hpp:100
ProSHADE_settings::setOutputFilename
void setOutputFilename(std::string oFileName)
Sets the requested output file name in the appropriate variable.
Definition: ProSHADE.cpp:636
ProSHADE_settings::useBiCubicInterpolationOnPeaks
bool useBiCubicInterpolationOnPeaks
This variable switch decides whether best symmetry is detected from peak indices, or whether bicubic ...
Definition: ProSHADE_settings.hpp:177
ProSHADE_run::getSymmetryAxis
std::vector< std::string > getSymmetryAxis(proshade_unsign axisNo)
This function returns a single symmetry axis as a vector of strings from the recommended symmetry axe...
Definition: ProSHADE.cpp:2361
ProSHADE_settings::setAxisComparisonThreshold
void setAxisComparisonThreshold(proshade_double axThres)
Sets the threshold for matching symmetry axes.
Definition: ProSHADE.cpp:955
ProSHADE_settings::forceP1
bool forceP1
Should the P1 spacegroup be forced on the input PDB files?
Definition: ProSHADE_settings.hpp:93
ProSHADE_settings::setTypicalNoiseSize
void setTypicalNoiseSize(proshade_single typNoi)
Sets the requested "fake" half-map kernel in the appropriate variable.
Definition: ProSHADE.cpp:505
ProSHADE_run
This class provides the access point to the library.
Definition: ProSHADE.hpp:39
ProSHADE_settings::task
ProSHADE_Task task
This custom type variable determines which task to perfom (i.e. symmetry detection,...
Definition: ProSHADE_settings.hpp:89
ProSHADE_settings::setMaskFilename
void setMaskFilename(std::string mskFln)
Sets where the mask should be saved.
Definition: ProSHADE.cpp:553
ProSHADE_run::getSymmetryFold
proshade_unsign getSymmetryFold(void)
This is the main accessor function for the user to get to know what symmetry fold ProSHADE has detect...
Definition: ProSHADE.cpp:1493
ProSHADE_settings::reBoxMap
bool reBoxMap
This switch decides whether re-boxing is needed.
Definition: ProSHADE_settings.hpp:136
ProSHADE_settings::setSymmetryRotFunPeaks
void setSymmetryRotFunPeaks(bool rotFunPeaks)
Sets the symmetry detection algorithm type.
Definition: ProSHADE.cpp:1134
ProSHADE_settings::normaliseMap
bool normaliseMap
Should the map be normalised to mean 0 sd 1?
Definition: ProSHADE_settings.hpp:120
ProSHADE_settings::addExtraSpace
proshade_single addExtraSpace
If this value is non-zero, this many angstroms of empty space will be added to the internal map.
Definition: ProSHADE_settings.hpp:146
ProSHADE_settings::setRequestedFold
void setRequestedFold(proshade_unsign val)
Sets the user requested symmetry fold.
Definition: ProSHADE.cpp:1060
ProSHADE_settings
This class stores all the settings and is passed to the executive classes instead of a multitude of p...
Definition: ProSHADE_settings.hpp:86
ProSHADE_run::getSymmetryType
std::string getSymmetryType(void)
This is the main accessor function for the user to get to know what symmetry type ProSHADE has detect...
Definition: ProSHADE.cpp:1483
ProSHADE_settings::setMapReboxing
void setMapReboxing(bool reBx)
Sets whether re-boxing needs to be done in the appropriate variable.
Definition: ProSHADE.cpp:569
ProSHADE_settings::setRecommendedSymmetry
void setRecommendedSymmetry(std::string val)
Sets the ProSHADE detected symmetry type.
Definition: ProSHADE.cpp:1009
ProSHADE_settings::setMaxSymmetryFold
void setMaxSymmetryFold(proshade_unsign maxFold)
Sets the maximum symmetry fold (well, the maximum prime symmetry fold).
Definition: ProSHADE.cpp:1162
ProSHADE_settings::correlationKernel
proshade_single correlationKernel
This value in Angstrom will be used as the kernel for the map-FHM correlation computation.
Definition: ProSHADE_settings.hpp:131
ProSHADE_settings::invertMap
bool invertMap
Should the map be inverted? Only use this if you think you have the wrong hand in your map.
Definition: ProSHADE_settings.hpp:123
ProSHADE_settings::setTaylorSeriesCap
void setTaylorSeriesCap(proshade_unsign tayCap)
Sets the requested Taylor series cap for the Gauss-Legendre integration in the appropriate variable.
Definition: ProSHADE.cpp:781
ProSHADE_settings::boundsSimilarityThreshold
proshade_signed boundsSimilarityThreshold
Number of indices which can be added just to make sure same size in indices is achieved.
Definition: ProSHADE_settings.hpp:138
ProSHADE_settings::setCorrelationMasking
void setCorrelationMasking(bool corMask)
Sets the requested map masking type in the appropriate variable.
Definition: ProSHADE.cpp:487
ProSHADE_settings::setProgressiveSphereMapping
void setProgressiveSphereMapping(bool progSphMap)
Sets the requested sphere mapping value settings approach in the appropriate variable.
Definition: ProSHADE.cpp:716
ProSHADE_settings::changeMapResolution
bool changeMapResolution
Should maps be re-sampled to obtain the required resolution?
Definition: ProSHADE_settings.hpp:99
ProSHADE_settings::setIntegrationOrder
void setIntegrationOrder(proshade_unsign intOrd)
Sets the requested order for the Gauss-Legendre integration in the appropriate variable.
Definition: ProSHADE.cpp:764
ProSHADE_settings::setRotationFunctionComputation
void setRotationFunctionComputation(bool rotfVal)
Sets whether the rotation function distance descriptor should be computed.
Definition: ProSHADE.cpp:832
ProSHADE_settings::rotationUncertainty
proshade_double rotationUncertainty
Alternative to bandwidth - the angle in degrees to which the rotation function accuracy should be com...
Definition: ProSHADE_settings.hpp:107
ProSHADE_settings::setOverlaySaveFile
void setOverlaySaveFile(std::string filename)
Sets the filename to which the overlay structure is to be save into.
Definition: ProSHADE.cpp:1106
ProSHADE_settings::overlayStructureName
std::string overlayStructureName
The filename to which the rotated and translated moving structure is to be saved.
Definition: ProSHADE_settings.hpp:181
ProSHADE_settings::moveToCOM
bool moveToCOM
Logical value stating whether the structure should be moved to have its Centre Of Mass (COM) in the m...
Definition: ProSHADE_settings.hpp:143
ProSHADE_settings::peakNeighbours
proshade_unsign peakNeighbours
Number of points in any direction that have to be lower than the considered index in order to conside...
Definition: ProSHADE_settings.hpp:161
ProSHADE_settings::smoothingFactor
proshade_double smoothingFactor
This factor decides how small the group sizes should be - larger factor means more smaller groups.
Definition: ProSHADE_settings.hpp:165
ProSHADE_settings::halfMapKernel
proshade_single halfMapKernel
This value in Angstrom will be used as the kernel for the "fake half-map" computation.
Definition: ProSHADE_settings.hpp:130
ProSHADE_internal_misc::checkMemoryAllocation
void checkMemoryAllocation(chVar checkVar, std::string fileP, unsigned int lineP, std::string funcP, std::string infoP="This error may occurs when ProSHADE requests memory to be\n : allocated to it and this operation fails. This could\n : happen when not enough memory is available, either due to\n : other processes using a lot of memory, or when the machine\n : does not have sufficient memory available. Re-run to see\n : if this problem persists.")
Checks if memory was allocated properly.
Definition: ProSHADE_misc.hpp:65
ProSHADE_settings::setGroupingSmoothingFactor
void setGroupingSmoothingFactor(proshade_double smFact)
Sets the grouping smoothing factor into the proper variable.
Definition: ProSHADE.cpp:921
ProSHADE_settings::setRequestedSymmetry
void setRequestedSymmetry(std::string val)
Sets the user requested symmetry type.
Definition: ProSHADE.cpp:1044
ProSHADE_settings::setExtraSpace
void setExtraSpace(proshade_single exSpace)
Sets the requested map extra space value in the appropriate variable.
Definition: ProSHADE.cpp:700
ProSHADE_settings::forceBounds
proshade_signed * forceBounds
These will be the boundaries to be forced upon the map.
Definition: ProSHADE_settings.hpp:140
ProSHADE_settings::recommendedSymmetryFold
proshade_unsign recommendedSymmetryFold
The fold of the recommended symmetry C or D type, 0 otherwise.
Definition: ProSHADE_settings.hpp:173
ProSHADE_settings::inputFiles
std::vector< std::string > inputFiles
This vector contains the filenames of all input structure files.
Definition: ProSHADE_settings.hpp:92
ProSHADE_settings::setMapResolutionChangeTriLinear
void setMapResolutionChangeTriLinear(bool mrChange)
Sets the requested map resolution change decision using tri-linear interpolation in the appropriate v...
Definition: ProSHADE.cpp:668
ProSHADE_settings::setMaskSaving
void setMaskSaving(bool savMsk)
Sets whether the mask should be saved.
Definition: ProSHADE.cpp:537
ProSHADE_settings::firstModelOnly
bool firstModelOnly
Shoud only the first PDB model be used, or should all models be used?
Definition: ProSHADE_settings.hpp:95
ProSHADE_settings::setBoundsThreshold
void setBoundsThreshold(proshade_signed boundsThres)
Sets the threshold for number of indices difference acceptable to make index sizes same in the approp...
Definition: ProSHADE.cpp:602
ProSHADE_settings::axisErrTolerance
proshade_double axisErrTolerance
Allowed error on vector axis in in dot product ( acos ( 1 - axErr ) is the allowed difference in radi...
Definition: ProSHADE_settings.hpp:169
ProSHADE_settings::useSameBounds
bool useSameBounds
Switch to say that the same boundaries as used for the first should be used for all input maps.
Definition: ProSHADE_settings.hpp:139
ProSHADE_settings::setMasking
void setMasking(bool mask)
Sets the requested map masking decision value in the appropriate variable.
Definition: ProSHADE.cpp:470
ProSHADE_settings::pdbBFactorNewVal
proshade_double pdbBFactorNewVal
Change all PDB B-factors to this value (for smooth maps).
Definition: ProSHADE_settings.hpp:103
ProSHADE_settings::setRecommendedFold
void setRecommendedFold(proshade_unsign val)
Sets the ProSHADE detected symmetry fold.
Definition: ProSHADE.cpp:1028
ProSHADE_settings::usePeakSearchInRotationFunctionSpace
bool usePeakSearchInRotationFunctionSpace
This variable switch decides whether symmetry detection will be done using peak search in rotation fu...
Definition: ProSHADE_settings.hpp:176
ProSHADE_settings::removeWaters
bool removeWaters
Should all waters be removed from input PDB files?
Definition: ProSHADE_settings.hpp:94
ProSHADE_settings::setResolution
void setResolution(proshade_single resolution)
Sets the requested resolution in the appropriate variable.
Definition: ProSHADE.cpp:357
ProSHADE_settings::printSettings
void printSettings(void)
This function prints the current values in the settings object.
Definition: ProSHADE.cpp:2098
ProSHADE_run::getNoStructures
proshade_unsign getNoStructures(void)
This function returns the number of structures used.
Definition: ProSHADE.cpp:2320
ProSHADE_settings::enLevMatrixPowerWeight
proshade_double enLevMatrixPowerWeight
If RRP matrices shell position is to be weighted by putting the position as an exponent,...
Definition: ProSHADE_settings.hpp:156
ProSHADE_settings::setNormalisation
void setNormalisation(bool normalise)
Sets the requested map normalisation value in the appropriate variable.
Definition: ProSHADE.cpp:389
ProSHADE_settings::setEnergyLevelsComputation
void setEnergyLevelsComputation(bool enLevDesc)
Sets whether the energy level distance descriptor should be computed.
Definition: ProSHADE.cpp:798
ProSHADE_settings::computeEnergyLevelsDesc
bool computeEnergyLevelsDesc
If true, the energy levels descriptor will be computed, otherwise all its computations will be omitte...
Definition: ProSHADE_settings.hpp:155