ProSHADE  0.7.6.1 (AUG 2021)
Protein Shape Detection
pyProSHADE.cpp
Go to the documentation of this file.
1 
22 //==================================================== Include PyBind11 header
23 #include <pybind11/pybind11.h>
24 #include <pybind11/stl.h>
25 #include <pybind11/numpy.h>
26 
27 //==================================================== Add the ProSHADE_settings and ProSHADE_run classes to the PyBind11 module
28 void add_settingsClass ( pybind11::module& pyProSHADE )
29 {
30  //================================================ Export the ProSHADE_settings class
31  pybind11::class_ < ProSHADE_settings > ( pyProSHADE, "ProSHADE_settings" )
32 
33  //============================================ Constructors (destructors do not need wrappers???)
34  .def ( pybind11::init < > ( ) )
35  .def ( pybind11::init < ProSHADE_Task > ( ), pybind11::arg ( "task" ) )
36 
37  //============================================ Member variables
38  .def_readwrite ( "task", &ProSHADE_settings::task )
39 
40  .def_readwrite ( "inputFiles", &ProSHADE_settings::inputFiles )
41  .def_readwrite ( "forceP1", &ProSHADE_settings::forceP1 )
42  .def_readwrite ( "removeNegativeDensity", &ProSHADE_settings::removeNegativeDensity )
43  .def_readwrite ( "removeWaters", &ProSHADE_settings::removeWaters )
44  .def_readwrite ( "firstModelOnly", &ProSHADE_settings::firstModelOnly )
45 
46  .def_readwrite ( "requestedResolution", &ProSHADE_settings::requestedResolution )
47  .def_readwrite ( "changeMapResolution", &ProSHADE_settings::changeMapResolution )
48  .def_readwrite ( "changeMapResolutionTriLinear", &ProSHADE_settings::changeMapResolutionTriLinear )
49 
50  .def_readwrite ( "pdbBFactorNewVal", &ProSHADE_settings::pdbBFactorNewVal )
51 
52  .def_readwrite ( "maxBandwidth", &ProSHADE_settings::maxBandwidth )
53  .def_readwrite ( "rotationUncertainty", &ProSHADE_settings::rotationUncertainty )
54 
55  .def_readwrite ( "usePhase", &ProSHADE_settings::usePhase )
56 
57  .def_readwrite ( "maxSphereDists", &ProSHADE_settings::maxSphereDists )
58 
59  .def_readwrite ( "integOrder", &ProSHADE_settings::integOrder )
60  .def_readwrite ( "taylorSeriesCap", &ProSHADE_settings::taylorSeriesCap )
61 
62  .def_readwrite ( "normaliseMap", &ProSHADE_settings::normaliseMap )
63 
64  .def_readwrite ( "invertMap", &ProSHADE_settings::invertMap )
65 
66  .def_readwrite ( "blurFactor", &ProSHADE_settings::blurFactor )
67  .def_readwrite ( "maskingThresholdIQRs", &ProSHADE_settings::maskingThresholdIQRs )
68  .def_readwrite ( "maskMap", &ProSHADE_settings::maskMap )
69  .def_readwrite ( "useCorrelationMasking", &ProSHADE_settings::useCorrelationMasking )
70  .def_readwrite ( "halfMapKernel", &ProSHADE_settings::halfMapKernel )
71  .def_readwrite ( "correlationKernel", &ProSHADE_settings::correlationKernel )
72  .def_readwrite ( "saveMask", &ProSHADE_settings::saveMask )
73  .def_readwrite ( "maskFileName", &ProSHADE_settings::maskFileName )
74  .def_readwrite ( "appliedMaskFileName", &ProSHADE_settings::appliedMaskFileName )
75 
76  .def_readwrite ( "fourierWeightsFileName", &ProSHADE_settings::fourierWeightsFileName )
77 
78  .def_readwrite ( "reBoxMap", &ProSHADE_settings::reBoxMap )
79  .def_readwrite ( "boundsExtraSpace", &ProSHADE_settings::boundsExtraSpace )
80  .def_readwrite ( "boundsSimilarityThreshold", &ProSHADE_settings::boundsSimilarityThreshold )
81  .def_readwrite ( "useSameBounds", &ProSHADE_settings::useSameBounds )
82  .def_readwrite ( "forceBounds", &ProSHADE_settings::forceBounds )
83 
84  .def_readwrite ( "moveToCOM", &ProSHADE_settings::moveToCOM )
85 
86  .def_readwrite ( "addExtraSpace", &ProSHADE_settings::addExtraSpace )
87 
88  .def_readwrite ( "progressiveSphereMapping", &ProSHADE_settings::progressiveSphereMapping )
89 
90  .def_readwrite ( "outName", &ProSHADE_settings::outName )
91 
92  .def_readwrite ( "computeEnergyLevelsDesc", &ProSHADE_settings::computeEnergyLevelsDesc )
93  .def_readwrite ( "enLevMatrixPowerWeight", &ProSHADE_settings::enLevMatrixPowerWeight )
94  .def_readwrite ( "computeTraceSigmaDesc", &ProSHADE_settings::computeTraceSigmaDesc )
95  .def_readwrite ( "computeRotationFuncDesc", &ProSHADE_settings::computeRotationFuncDesc )
96 
97  .def_readwrite ( "peakNeighbours", &ProSHADE_settings::peakNeighbours )
98  .def_readwrite ( "noIQRsFromMedianNaivePeak", &ProSHADE_settings::noIQRsFromMedianNaivePeak )
99 
100  .def_readwrite ( "smoothingFactor", &ProSHADE_settings::smoothingFactor )
101 
102  .def_readwrite ( "symMissPeakThres", &ProSHADE_settings::symMissPeakThres )
103  .def_readwrite ( "axisErrTolerance", &ProSHADE_settings::axisErrTolerance )
104  .def_readwrite ( "axisErrToleranceDefault", &ProSHADE_settings::axisErrToleranceDefault )
105  .def_readwrite ( "minSymPeak", &ProSHADE_settings::minSymPeak )
106  .def_readwrite ( "recommendedSymmetryType", &ProSHADE_settings::recommendedSymmetryType )
107  .def_readwrite ( "recommendedSymmetryFold", &ProSHADE_settings::recommendedSymmetryFold )
108  .def_readwrite ( "requestedSymmetryType", &ProSHADE_settings::requestedSymmetryType )
109  .def_readwrite ( "requestedSymmetryFold", &ProSHADE_settings::requestedSymmetryFold )
110  .def_readwrite ( "useBiCubicInterpolationOnPeaks", &ProSHADE_settings::useBiCubicInterpolationOnPeaks )
111  .def_readwrite ( "maxSymmetryFold", &ProSHADE_settings::maxSymmetryFold )
112  .def_readwrite ( "fscThreshold", &ProSHADE_settings::fscThreshold )
113  .def_readwrite ( "peakThresholdMin", &ProSHADE_settings::peakThresholdMin )
114 
115  .def_readwrite ( "overlayStructureName", &ProSHADE_settings::overlayStructureName )
116  .def_readwrite ( "rotTrsJSONFile", &ProSHADE_settings::rotTrsJSONFile )
117 
118  .def_readwrite ( "verbose", &ProSHADE_settings::verbose )
119 
120  //============================================ Mutators
121  .def ( "addStructure", &ProSHADE_settings::addStructure, "Adds a structure file name to the appropriate variable.", pybind11::arg ( "structure" ) )
122  .def ( "setResolution", &ProSHADE_settings::setResolution, "This function sets the resolution in the appropriate variable.", pybind11::arg ( "resolution" ) )
123  .def ( "setPDBBFactor", &ProSHADE_settings::setPDBBFactor, "Sets the requested B-factor value for PDB files in the appropriate variable.", pybind11::arg ( "newBF" ) )
124  .def ( "setNormalisation", &ProSHADE_settings::setNormalisation, "Sets the requested map normalisation value in the appropriate variable.", pybind11::arg ( "normalise" ) )
125  .def ( "setMapInversion", &ProSHADE_settings::setMapInversion, "Sets the requested map inversion value in the appropriate variable.", pybind11::arg ( "mInv" ) )
126  .def ( "setVerbosity", &ProSHADE_settings::setVerbosity, "Sets the requested verbosity in the appropriate variable.", pybind11::arg ( "verbosity" ) )
127  .def ( "setMaskBlurFactor", &ProSHADE_settings::setMaskBlurFactor, "Sets the requested map blurring factor in the appropriate variable.", pybind11::arg ( "blurFac" ) )
128  .def ( "setMaskIQR", &ProSHADE_settings::setMaskIQR, "Sets the requested number of IQRs for masking threshold in the appropriate variable.", pybind11::arg ( "noIQRs" ) )
129  .def ( "setMasking", &ProSHADE_settings::setMasking, "Sets the requested map masking decision value in the appropriate variable.", pybind11::arg ( "mask" ) )
130  .def ( "setCorrelationMasking", &ProSHADE_settings::setCorrelationMasking, "Sets the requested map masking type in the appropriate variable.", pybind11::arg ( "corMask" ) )
131  .def ( "setTypicalNoiseSize", &ProSHADE_settings::setTypicalNoiseSize, "Sets the requested \"fake\" half-map kernel in the appropriate variable.", pybind11::arg ( "typNoi" ) )
132  .def ( "setMinimumMaskSize", &ProSHADE_settings::setMinimumMaskSize, "Sets the requested minimum mask size.", pybind11::arg ( "minMS" ) )
133  .def ( "setMaskSaving", &ProSHADE_settings::setMaskSaving, "Sets whether the mask should be saved.", pybind11::arg ( "savMsk" ) )
134  .def ( "setMaskFilename", &ProSHADE_settings::setMaskFilename, "Sets where the mask should be saved.", pybind11::arg ( "mskFln" ) )
135  .def ( "setAppliedMaskFilename", &ProSHADE_settings::setAppliedMaskFilename, "Sets the filename of the mask data that should be applied to the input map.", pybind11::arg ( "mskFln" ) )
136  .def ( "setFourierWeightsFilename", &ProSHADE_settings::setFourierWeightsFilename, "Sets the filename of the Fourier weights data that should be applied to the input map.", pybind11::arg ( "fWgFln" ) )
137  .def ( "setMapReboxing", &ProSHADE_settings::setMapReboxing, "Sets whether re-boxing needs to be done in the appropriate variable.", pybind11::arg ( "reBx" ) )
138  .def ( "setBoundsSpace", &ProSHADE_settings::setBoundsSpace, "Sets the requested number of angstroms for extra space in re-boxing in the appropriate variable.", pybind11::arg ( "boundsExSp" ) )
139  .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" ) )
140  .def ( "setSameBoundaries", &ProSHADE_settings::setSameBoundaries, "Sets whether same boundaries should be used in the appropriate variable.", pybind11::arg ( "sameB" ) )
141  .def ( "setOutputFilename", &ProSHADE_settings::setOutputFilename, "Sets the requested output file name in the appropriate variable.", pybind11::arg ( "oFileName" ) )
142  .def ( "setMapResolutionChange", &ProSHADE_settings::setMapResolutionChange, "Sets the requested map resolution change decision in the appropriate variable.", pybind11::arg ( "mrChange" ) )
143  .def ( "setMapResolutionChangeTriLinear", &ProSHADE_settings::setMapResolutionChangeTriLinear, "Sets the requested map resolution change decision using tri-linear interpolation in the appropriate variable.", pybind11::arg ( "mrChange" ) )
144  .def ( "setMapCentering", &ProSHADE_settings::setMapCentering, "Sets the requested map centering decision value in the appropriate variable.", pybind11::arg ( "com" ) )
145  .def ( "setExtraSpace", &ProSHADE_settings::setExtraSpace, "Sets the requested map extra space value in the appropriate variable.", pybind11::arg ( "exSpace" ) )
146  .def ( "setBandwidth", &ProSHADE_settings::setBandwidth, "Sets the requested spherical harmonics bandwidth in the appropriate variable.", pybind11::arg ( "band" ) )
147  .def ( "setProgressiveSphereMapping", &ProSHADE_settings::setProgressiveSphereMapping, "Sets the requested sphere mapping value settings approach in the appropriate variable.", pybind11::arg ( "progSphMap" ) )
148  .def ( "setSphereDistances", &ProSHADE_settings::setSphereDistances, "Sets the requested distance between spheres in the appropriate variable.", pybind11::arg ( "sphDist" ) )
149  .def ( "setIntegrationOrder", &ProSHADE_settings::setIntegrationOrder, "Sets the requested order for the Gauss-Legendre integration in the appropriate variable.", pybind11::arg ( "intOrd" ) )
150  .def ( "setTaylorSeriesCap", &ProSHADE_settings::setTaylorSeriesCap, "Sets the requested Taylor series cap for the Gauss-Legendre integration in the appropriate variable.", pybind11::arg ( "tayCap" ) )
151  .def ( "setEnergyLevelsComputation", &ProSHADE_settings::setEnergyLevelsComputation, "Sets whether the energy level distance descriptor should be computed.", pybind11::arg ( "enLevDesc" ) )
152  .def ( "setTraceSigmaComputation", &ProSHADE_settings::setTraceSigmaComputation, "Sets whether the trace sigma distance descriptor should be computed.", pybind11::arg ( "trSigVal" ) )
153  .def ( "setRotationFunctionComputation", &ProSHADE_settings::setRotationFunctionComputation, "Sets whether the rotation function distance descriptor should be computed.", pybind11::arg ( "rotfVal" ) )
154  .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" ) )
155  .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" ) )
156  .def ( "setPhaseUsage", &ProSHADE_settings::setPhaseUsage, "Sets whether the phase information will be used.", pybind11::arg ( "phaseUsage" ) )
157  .def ( "setEnLevShellWeight", &ProSHADE_settings::setEnLevShellWeight, "Sets the weight of shell position for the energy levels computation.", pybind11::arg ( "mPower" ) )
158  .def ( "setGroupingSmoothingFactor", &ProSHADE_settings::setGroupingSmoothingFactor, "Sets the grouping smoothing factor into the proper variable.", pybind11::arg ( "smFact" ) )
159  .def ( "setMissingPeakThreshold", &ProSHADE_settings::setMissingPeakThreshold, "Sets the threshold for starting the missing peaks procedure.", pybind11::arg ( "mpThres" ) )
160  .def ( "setAxisComparisonThreshold", &ProSHADE_settings::setAxisComparisonThreshold, "Sets the threshold for matching symmetry axes.", pybind11::arg ( "axThres" ) )
161  .def ( "setAxisComparisonThresholdBehaviour", &ProSHADE_settings::setAxisComparisonThresholdBehaviour, "Sets the automatic symmetry axis tolerance decreasing.", pybind11::arg ( "behav" ) )
162  .def ( "setMinimumPeakForAxis", &ProSHADE_settings::setMinimumPeakForAxis, "Sets the minimum peak height for symmetry axis to be considered.", pybind11::arg ( "minSP" ) )
163  .def ( "setRecommendedSymmetry", &ProSHADE_settings::setRecommendedSymmetry, "Sets the ProSHADE detected symmetry type.", pybind11::arg ( "val" ) )
164  .def ( "setRecommendedFold", &ProSHADE_settings::setRecommendedFold, "Sets the ProSHADE detected symmetry fold.", pybind11::arg ( "val" ) )
165  .def ( "setRequestedSymmetry", &ProSHADE_settings::setRequestedSymmetry, "Sets the user requested symmetry type.", pybind11::arg ( "val" ) )
166  .def ( "setRequestedFold", &ProSHADE_settings::setRequestedFold, "Sets the user requested symmetry fold.", pybind11::arg ( "val" ) )
167  .def ( "setDetectedSymmetry", &ProSHADE_settings::setDetectedSymmetry, "Sets the final detected symmetry axes information.", pybind11::arg ( "sym" ) )
168  .def ( "setOverlaySaveFile", &ProSHADE_settings::setOverlaySaveFile, "Sets the filename to which the overlay structure is to be save into.", pybind11::arg ( "filename" ) )
169  .def ( "setOverlayJsonFile", &ProSHADE_settings::setOverlayJsonFile, "Sets the filename to which the overlay operations are to be save into.", pybind11::arg ( "filename" ) )
170  .def ( "setBicubicInterpolationSearch", &ProSHADE_settings::setBicubicInterpolationSearch, "Sets the bicubic interpolation on peaks.", pybind11::arg ( "bicubPeaks" ) )
171  .def ( "setMaxSymmetryFold", &ProSHADE_settings::setMaxSymmetryFold, "Sets the maximum symmetry fold (well, the maximum prime symmetry fold).", pybind11::arg ( "maxFold" ) )
172  .def ( "setFSCThreshold", &ProSHADE_settings::setFSCThreshold, "Sets the minimum FSC threshold for axis to be considered detected.", pybind11::arg ( "fscThr" ) )
173  .def ( "setPeakThreshold", &ProSHADE_settings::setPeakThreshold, "Sets the minimum peak height threshold for axis to be considered possible.", pybind11::arg ( "peakThr" ) )
174  .def ( "setNegativeDensity", &ProSHADE_settings::setNegativeDensity, "Sets the internal variable deciding whether input files negative density should be removed.", pybind11::arg ( "nDens" ) )
175 
176  //============================================ Command line parsing
177  .def ( "getCommandLineParams",
178  [] ( ProSHADE_settings &self, std::vector < std::string > args )
179  {
180  std::vector < char * > cstrs; cstrs.reserve ( args.size() );
181 
182  for ( auto &s : args )
183  cstrs.push_back ( const_cast < char * > ( s.c_str ( ) ) );
184 
185  return self.getCommandLineParams ( static_cast< int > ( cstrs.size ( ) ), cstrs.data ( ) );
186  }, "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." )
187 
188  //============================================ Debugging
189  .def ( "printSettings", &ProSHADE_settings::printSettings, "This function prints the current values in the settings object." )
190 
191  //============================================ Description
192  .def ( "__repr__", [] ( ) { return "<ProSHADE_settings class object> (Settings class is used to set all settings values in a single place)"; } );
193 
194  //================================================ Export the ProSHADE_run class
195  pybind11::class_ < ProSHADE_run > ( pyProSHADE, "ProSHADE_run" )
196 
197  //============================================ Constructors (destructors do not need wrappers???)
198  .def ( pybind11::init < ProSHADE_settings* > ( ) )
199 
200  //============================================ General accessors
201  .def ( "getNoStructures", &ProSHADE_run::getNoStructures, "This function returns the number of structures used." )
202  .def ( "getVerbose", &ProSHADE_run::getVerbose, "This function returns the verbose value." )
203 
204  //============================================ Distances results accessor functions wrapped as lambda functions for numpy return types
205  .def ( "getEnergyLevelsVector",
206  [] ( ProSHADE_run &self ) -> pybind11::array_t < float >
207  {
208  //== Get the values
209  std::vector< proshade_double > vals = self.getEnergyLevelsVector ();
210 
211  //== Allocate memory for the numpy values
212  float* npVals = new float[static_cast<unsigned int> (vals.size())];
213  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
214 
215  //== Copy values
216  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = static_cast< float > ( vals.at(iter) ); }
217 
218  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
219  pybind11::capsule pyCapsuleEnLevs ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
220 
221  //== Copy the value
222  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<unsigned int> (vals.size()) }, // Shape
223  { sizeof(float) }, // C-stype strides
224  npVals, // Data
225  pyCapsuleEnLevs ); // Capsule (C++ destructor, basically)
226 
227  //== Done
228  return ( retArr );
229  }, "This function returns the energy level distances vector from the first to all other structures." )
230 
231  .def ( "getTraceSigmaVector",
232  [] ( ProSHADE_run &self ) -> pybind11::array_t < float >
233  {
234  //== Get the values
235  std::vector< proshade_double > vals = self.getTraceSigmaVector ();
236 
237  //== Allocate memory for the numpy values
238  float* npVals = new float[static_cast<unsigned int> (vals.size())];
239  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
240 
241  //== Copy values
242  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = static_cast< float > ( vals.at(iter) ); }
243 
244  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
245  pybind11::capsule pyCapsuleTrSigs ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
246 
247  //== Copy the value
248  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<unsigned int> (vals.size()) }, // Shape
249  { sizeof(float) }, // C-stype strides
250  npVals, // Data
251  pyCapsuleTrSigs ); // Capsule (C++ destructor, basically)
252 
253  //== Done
254  return ( retArr );
255  }, "This function returns the trace sigma distances vector from the first to all other structures." )
256 
257  .def ( "getRotationFunctionVector",
258  [] ( ProSHADE_run &self ) -> pybind11::array_t < float >
259  {
260  //== Get the values
261  std::vector< proshade_double > vals = self.getRotationFunctionVector ();
262 
263  //== Allocate memory for the numpy values
264  float* npVals = new float[static_cast<unsigned int> (vals.size())];
265  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
266 
267  //== Copy values
268  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = static_cast< float > ( vals.at(iter) ); }
269 
270  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
271  pybind11::capsule pyCapsuleRotFun ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
272 
273  //== Copy the value
274  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<unsigned int> (vals.size()) }, // Shape
275  { sizeof(float) }, // C-stype strides
276  npVals, // Data
277  pyCapsuleRotFun ); // Capsule (C++ destructor, basically)
278 
279  //== Done
280  return ( retArr );
281  }, "This function returns the full rotation function distances vector from the first to all other structures." )
282 
283  //============================================ Symmetry results accessor functions
284  .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." )
285  .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." )
286  .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" ) )
287  .def ( "getAllCSyms",
288  [] ( ProSHADE_run &self ) -> pybind11::array_t < float >
289  {
290  //== Get the values
291  std::vector< std::vector< proshade_double > > vals = self.getAllCSyms ();
292 
293  //== Allocate memory for the numpy values
294  float* npVals = new float[static_cast<unsigned int> ( vals.size() * 7 )];
295  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
296 
297  //== Copy values
298  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { for ( proshade_unsign it = 0; it < 7; it++ ) { npVals[(iter*7)+it] = static_cast< float > ( vals.at(iter).at(it) ); } }
299 
300  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
301  pybind11::capsule pyCapsuleSymList ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
302 
303  //== Copy the value
304  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<int> ( vals.size() ), static_cast<int> ( 7 ) }, // Shape
305  { 7 * sizeof(float), sizeof(float) }, // C-stype strides
306  npVals, // Data
307  pyCapsuleSymList ); // Capsule
308 
309  //== Done
310  return ( retArr );
311  }, "This function returns a all symmetry axes as a 2D numpy array." )
312  .def ( "getMapCOMProcessChange",
313  [] ( ProSHADE_run &self ) -> pybind11::array_t < float >
314  {
315  //== Get the values
316  std::vector< proshade_double > vals = self.getMapCOMProcessChange ();
317 
318  //== Allocate memory for the numpy values
319  float* npVals = new float[static_cast<unsigned int> ( 3 )];
320  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
321 
322  //== Copy values
323  for ( proshade_unsign iter = 0; iter < 3; iter++ ) { npVals[iter] = static_cast< float > ( vals.at(iter) ); }
324 
325  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
326  pybind11::capsule pyCapsuleSymShift ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
327 
328  //== Copy the value
329  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<int> ( vals.size() ) }, // Shape
330  { sizeof(float) }, // C-stype strides
331  npVals, // Data
332  pyCapsuleSymShift ); // Capsule
333 
334  //== Done
335  return ( retArr );
336  }, "This function returns the shift in Angstrom applied to the internal map representation in order to align its COM with the centre of box." )
337 
338  //============================================ Reboxing results accessor functions as lambda functions directly returning numpy arrays
339  .def ( "getOriginalBounds",
340  [] ( ProSHADE_run &self, proshade_unsign strNo ) -> pybind11::array_t < float >
341  {
342  //== Get values
343  std::vector< proshade_signed > vals = self.getOriginalBounds ( strNo );
344 
345  //== Allocate memory for the numpy values
346  float* npVals = new float[static_cast<proshade_unsign> ( vals.size() )];
347  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
348 
349  //== Copy values
350  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = vals.at(iter); }
351 
352  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
353  pybind11::capsule pyCapsuleOrigBnds ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
354 
355  //== Copy the value
356  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<proshade_unsign> ( vals.size() ) }, // Shape
357  { sizeof(float) }, // C-stype strides
358  npVals, // Data
359  pyCapsuleOrigBnds ); // Capsule (C++ destructor, basically)
360 
361  //== Done
362  return ( retArr );
363  }, "This function returns the original structure boundaries as numpy array." )
364 
365  .def ( "getReBoxedBounds",
366  [] ( ProSHADE_run &self, proshade_unsign strNo ) -> pybind11::array_t < float >
367  {
368  //== Get values
369  std::vector< proshade_signed > vals = self.getReBoxedBounds ( strNo );
370 
371  //== Allocate memory for the numpy values
372  float* npVals = new float[static_cast<proshade_unsign> ( vals.size() )];
373  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
374 
375  //== Copy values
376  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = vals.at(iter); }
377 
378  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
379  pybind11::capsule pyCapsuleReBoBnds ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
380 
381  //== Copy the value
382  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<proshade_unsign> ( vals.size() ) }, // Shape
383  { sizeof(float) }, // C-stype strides
384  npVals, // Data
385  pyCapsuleReBoBnds ); // Capsule (C++ destructor, basically)
386 
387  //== Done
388  return ( retArr );
389  }, "This function returns the re-boxed structure boundaries as numpy array." )
390 
391 
392  .def ( "getReBoxedMap",
393  [] ( ProSHADE_run &self, proshade_unsign strNo ) -> pybind11::array_t < float >
394  {
395  //== Get the values
396  std::vector< proshade_signed > vals = self.getReBoxedBounds ( strNo );
397 
398  //== Determine dimensions
399  proshade_unsign xDim = static_cast< proshade_unsign > ( vals.at(1) ) - static_cast< proshade_unsign > ( vals.at(0) ) + 1;
400  proshade_unsign yDim = static_cast< proshade_unsign > ( vals.at(3) ) - static_cast< proshade_unsign > ( vals.at(2) ) + 1;
401  proshade_unsign zDim = static_cast< proshade_unsign > ( vals.at(5) ) - static_cast< proshade_unsign > ( vals.at(4) ) + 1;
402 
403  //== Allocate memory for the numpy values
404  float* npVals = new float[xDim * yDim * zDim];
405  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
406 
407  //== Copy values
408  for ( proshade_unsign iter = 0; iter < (xDim * yDim * zDim); iter++ ) { npVals[iter] = static_cast< float > ( self.getMapValue ( strNo, iter ) ); }
409 
410  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
411  pybind11::capsule pyCapsuleRebMap ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
412 
413  //== Copy the value
414  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { xDim, yDim, zDim }, // Shape
415  { yDim * zDim * sizeof(float), zDim * sizeof(float), sizeof(float) }, // C-stype strides
416  npVals, // Data
417  pyCapsuleRebMap ); // Capsule
418 
419  //== Done
420  return ( retArr );
421  }, "This function returns the re-boxed structure map as a numpy 3D array." )
422 
423  //============================================ Overlay results accessor functions as lambda functions directly returning numpy arrays
424  .def ( "getEulerAngles",
425  [] ( ProSHADE_run &self ) -> pybind11::array_t < float >
426  {
427  //== Get the values
428  std::vector< proshade_double > vals = self.getEulerAngles ( );
429 
430  //== Allocate memory for the numpy values
431  float* npVals = new float[static_cast<proshade_unsign> ( vals.size() )];
432  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
433 
434  //== Copy values
435  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = static_cast< float > ( vals.at(iter) ); }
436 
437  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
438  pybind11::capsule pyCapsuleEulAngs ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
439 
440  //== Copy the value
441  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<proshade_unsign> ( vals.size() ) }, // Shape
442  { sizeof(float) }, // C-stype strides
443  npVals, // Data
444  pyCapsuleEulAngs ); // Capsule
445 
446  //== Done
447  return ( retArr );
448  }, "This function returns the vector of Euler angles with best overlay correlation." )
449 
450  .def ( "getOptimalRotMat",
451  [] ( ProSHADE_run &self ) -> pybind11::array_t < float >
452  {
453  //== Get the values
454  std::vector< proshade_double > vals = self.getOptimalRotMat ( );
455 
456  //== Allocate memory for the numpy values
457  float* npVals = new float[static_cast<proshade_unsign> ( vals.size() )];
458  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
459 
460  //== Copy values
461  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = static_cast< float > ( vals.at(iter) ); }
462 
463  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
464  pybind11::capsule pyCapsuleRotMat ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
465 
466  //== Copy the value
467  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { 3, 3 }, // Shape
468  { 3 * sizeof(float), sizeof(float) }, // C-stype strides
469  npVals, // Data
470  pyCapsuleRotMat ); // Capsule
471 
472  //== Done
473  return ( retArr );
474  }, "This function returns the vector of Euler angles with best overlay correlation." )
475 
476  .def ( "getTranslationToOrigin",
477  [] ( ProSHADE_run &self ) -> pybind11::array_t < float >
478  {
479  //== Get the values
480  std::vector< proshade_double > vals = self.getTranslationToOrigin ( );
481 
482  //== Allocate memory for the numpy values
483  float* npVals = new float[static_cast<proshade_unsign> ( vals.size() )];
484  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
485 
486  //== Copy values
487  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = static_cast< float > ( vals.at(iter) ); }
488 
489  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
490  pybind11::capsule pyCapsuleTTO ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
491 
492  //== Copy the value
493  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<proshade_unsign> ( vals.size() ) }, // Shape
494  { sizeof(float) }, // C-stype strides
495  npVals, // Data
496  pyCapsuleTTO ); // Capsule
497 
498  //== Done
499  return ( retArr );
500  }, "This function returns the negative values of the position of the rotation centre (the point about which the rotation should be done)." )
501  .def ( "getOriginToOverlayTranslation",
502  [] ( ProSHADE_run &self ) -> pybind11::array_t < float >
503  {
504  //== Get the values
505  std::vector< proshade_double > vals = self.getOriginToOverlayTranslation ( );
506 
507  //== Allocate memory for the numpy values
508  float* npVals = new float[static_cast<proshade_unsign> ( vals.size() )];
509  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
510 
511  //== Copy values
512  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = static_cast< float > ( vals.at(iter) ); }
513 
514  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
515  pybind11::capsule pyCapsuleOTOT ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
516 
517  //== Copy the value
518  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<proshade_unsign> ( vals.size() ) }, // Shape
519  { sizeof(float) }, // C-stype strides
520  npVals, // Data
521  pyCapsuleOTOT ); // Capsule
522 
523  //== Done
524  return ( retArr );
525  }, "This function returns the translation required to move the structure from origin to optimal overlay." )
526 
527 
528  //============================================ Description
529  .def ( "__repr__", [] ( ) { return "<ProSHADE_run class object> (Run class constructor takes a ProSHADE_settings object and completes a single run according to the settings object information)"; } );
530 }
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:118
ProSHADE_settings::setOverlayJsonFile
void setOverlayJsonFile(std::string filename)
Sets the filename to which the overlay operations are to be save into.
Definition: ProSHADE.cpp:1367
ProSHADE_settings::maxBandwidth
proshade_unsign maxBandwidth
The bandwidth of spherical harmonics decomposition for the largest sphere.
Definition: ProSHADE_settings.hpp:58
ProSHADE_settings::integOrder
proshade_unsign integOrder
The order required for full Gauss-Legendre integration between the spheres.
Definition: ProSHADE_settings.hpp:68
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:128
ProSHADE_settings::rotTrsJSONFile
std::string rotTrsJSONFile
The filename to which the rotation and translation operations are to be saved into.
Definition: ProSHADE_settings.hpp:139
ProSHADE_settings::setEnLevShellWeight
void setEnLevShellWeight(proshade_double mPower)
Sets the weight of shell position for the energy levels computation.
Definition: ProSHADE.cpp:1101
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:113
ProSHADE_settings::setTraceSigmaComputation
void setTraceSigmaComputation(bool trSigVal)
Sets whether the trace sigma distance descriptor should be computed.
Definition: ProSHADE.cpp:993
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:114
ProSHADE_settings::setSameBoundaries
void setSameBoundaries(bool sameB)
Sets whether same boundaries should be used in the appropriate variable.
Definition: ProSHADE.cpp:749
ProSHADE_settings::maxSymmetryFold
proshade_unsign maxSymmetryFold
The highest symmetry fold to search for.
Definition: ProSHADE_settings.hpp:133
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:69
ProSHADE_settings::setMinimumMaskSize
void setMinimumMaskSize(proshade_single minMS)
Sets the requested minimum mask size.
Definition: ProSHADE.cpp:587
ProSHADE_settings::setAppliedMaskFilename
void setAppliedMaskFilename(std::string mskFln)
Sets the filename of the mask data that should be applied to the input map.
Definition: ProSHADE.cpp:647
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:93
ProSHADE_settings::setBicubicInterpolationSearch
void setBicubicInterpolationSearch(bool bicubPeaks)
Sets the bicubic interpolation on peaks.
Definition: ProSHADE.cpp:1385
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:1057
ProSHADE_settings::outName
std::string outName
The file name where the output structure(s) should be saved.
Definition: ProSHADE_settings.hpp:108
ProSHADE_settings::setMapInversion
void setMapInversion(bool mInv)
Sets the requested map inversion value in the appropriate variable.
Definition: ProSHADE.cpp:443
ProSHADE_settings::setMapResolutionChange
void setMapResolutionChange(bool mrChange)
Sets the requested map resolution change decision in the appropriate variable.
Definition: ProSHADE.cpp:790
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:1035
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:78
ProSHADE_settings::maskMap
bool maskMap
Should the map be masked from noise?
Definition: ProSHADE_settings.hpp:80
ProSHADE_settings::requestedSymmetryFold
proshade_unsign requestedSymmetryFold
The fold of the requested symmetry (only applicable to C and D symmetry types).
Definition: ProSHADE_settings.hpp:131
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:130
ProSHADE_settings::setSphereDistances
void setSphereDistances(proshade_single sphDist)
Sets the requested distance between spheres in the appropriate variable.
Definition: ProSHADE.cpp:890
ProSHADE_settings::setVerbosity
void setVerbosity(proshade_signed verbosity)
Sets the requested verbosity in the appropriate variable.
Definition: ProSHADE.cpp:463
ProSHADE_settings::removeNegativeDensity
bool removeNegativeDensity
Should the negative density be removed from input files?
Definition: ProSHADE_settings.hpp:47
ProSHADE_settings::setMinimumPeakForAxis
void setMinimumPeakForAxis(proshade_double minSP)
Sets the minimum peak height for symmetry axis to be considered.
Definition: ProSHADE.cpp:1209
ProSHADE_settings::saveMask
bool saveMask
Should the mask be saved?
Definition: ProSHADE_settings.hpp:84
ProSHADE_settings::requestedResolution
proshade_single requestedResolution
The resolution to which the calculations are to be done.
Definition: ProSHADE_settings.hpp:50
ProSHADE_settings::minSymPeak
proshade_double minSymPeak
Minimum average peak for symmetry axis to be considered as "real".
Definition: ProSHADE_settings.hpp:127
ProSHADE_settings::setAxisComparisonThresholdBehaviour
void setAxisComparisonThresholdBehaviour(bool behav)
Sets the automatic symmetry axis tolerance decreasing.
Definition: ProSHADE.cpp:1188
ProSHADE_settings::setDetectedSymmetry
void setDetectedSymmetry(proshade_double *sym)
Sets the final detected symmetry axes information.
Definition: ProSHADE.cpp:1315
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:105
ProSHADE_settings::maxSphereDists
proshade_single maxSphereDists
The distance between spheres in spherical mapping for the largest sphere.
Definition: ProSHADE_settings.hpp:65
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:124
ProSHADE_settings::setPDBBFactor
void setPDBBFactor(proshade_double newBF)
Sets the requested B-factor value for PDB files in the appropriate variable.
Definition: ProSHADE.cpp:403
ProSHADE_settings::addStructure
void addStructure(std::string structure)
Adds a structure file name to the appropriate variable.
Definition: ProSHADE.cpp:363
ProSHADE_settings::setMaskIQR
void setMaskIQR(proshade_single noIQRs)
Sets the requested number of IQRs for masking threshold in the appropriate variable.
Definition: ProSHADE.cpp:504
ProSHADE_settings::setMissingPeakThreshold
void setMissingPeakThreshold(proshade_double mpThres)
Sets the threshold for starting the missing peaks procedure.
Definition: ProSHADE.cpp:1144
ProSHADE_run::getVerbose
proshade_signed getVerbose(void)
This function returns the verbose value.
Definition: ProSHADE.cpp:2792
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:79
ProSHADE_settings::useCorrelationMasking
bool useCorrelationMasking
Should the blurring masking (false) or the correlation masking (true) be used?
Definition: ProSHADE_settings.hpp:81
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:62
ProSHADE_settings::setMapCentering
void setMapCentering(bool com)
Sets the requested map centering decision value in the appropriate variable.
Definition: ProSHADE.cpp:830
ProSHADE_settings::setMaskBlurFactor
void setMaskBlurFactor(proshade_single blurFac)
Sets the requested map blurring factor in the appropriate variable.
Definition: ProSHADE.cpp:483
ProSHADE_settings::setPeakThreshold
void setPeakThreshold(proshade_double peakThr)
Sets the minimum peak height threshold for axis to be considered possible.
Definition: ProSHADE.cpp:1439
ProSHADE_settings::setPhaseUsage
void setPhaseUsage(bool phaseUsage)
Sets whether the phase information will be used.
Definition: ProSHADE.cpp:1079
ProSHADE_settings::maskFileName
std::string maskFileName
The filename to which mask should be saved.
Definition: ProSHADE_settings.hpp:85
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:142
ProSHADE_settings::setFSCThreshold
void setFSCThreshold(proshade_double fscThr)
Sets the minimum FSC threshold for axis to be considered detected.
Definition: ProSHADE.cpp:1421
ProSHADE_settings::setBandwidth
void setBandwidth(proshade_unsign band)
Sets the requested spherical harmonics bandwidth in the appropriate variable.
Definition: ProSHADE.cpp:870
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:708
ProSHADE_settings::changeMapResolutionTriLinear
bool changeMapResolutionTriLinear
Should maps be re-sampled to obtain the required resolution?
Definition: ProSHADE_settings.hpp:52
ProSHADE_settings::setOutputFilename
void setOutputFilename(std::string oFileName)
Sets the requested output file name in the appropriate variable.
Definition: ProSHADE.cpp:770
ProSHADE_settings::useBiCubicInterpolationOnPeaks
bool useBiCubicInterpolationOnPeaks
This variable switch decides whether best symmetry is detected from peak indices, or whether bicubic ...
Definition: ProSHADE_settings.hpp:132
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:2826
ProSHADE_settings::setAxisComparisonThreshold
void setAxisComparisonThreshold(proshade_double axThres)
Sets the threshold for matching symmetry axes.
Definition: ProSHADE.cpp:1165
ProSHADE_settings::peakThresholdMin
proshade_double peakThresholdMin
The threshold for peak height above which axes are considered possible.
Definition: ProSHADE_settings.hpp:135
ProSHADE_settings::forceP1
bool forceP1
Should the P1 spacegroup be forced on the input PDB files?
Definition: ProSHADE_settings.hpp:44
ProSHADE_settings::setTypicalNoiseSize
void setTypicalNoiseSize(proshade_single typNoi)
Sets the requested "fake" half-map kernel in the appropriate variable.
Definition: ProSHADE.cpp:567
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:40
ProSHADE_settings::setMaskFilename
void setMaskFilename(std::string mskFln)
Sets where the mask should be saved.
Definition: ProSHADE.cpp:627
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:1814
ProSHADE_settings::reBoxMap
bool reBoxMap
This switch decides whether re-boxing is needed.
Definition: ProSHADE_settings.hpp:92
ProSHADE_settings::appliedMaskFileName
std::string appliedMaskFileName
The filename from which mask data will be read from.
Definition: ProSHADE_settings.hpp:86
ProSHADE_settings::normaliseMap
bool normaliseMap
Should the map be normalised to mean 0 sd 1?
Definition: ProSHADE_settings.hpp:72
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:102
ProSHADE_settings::setRequestedFold
void setRequestedFold(proshade_unsign val)
Sets the user requested symmetry fold.
Definition: ProSHADE.cpp:1294
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:37
ProSHADE_settings::fscThreshold
proshade_double fscThreshold
The threshold for FSC value under which the axis is considered to be likely noise.
Definition: ProSHADE_settings.hpp:134
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:1799
ProSHADE_settings::setMapReboxing
void setMapReboxing(bool reBx)
Sets whether re-boxing needs to be done in the appropriate variable.
Definition: ProSHADE.cpp:687
ProSHADE_settings::setRecommendedSymmetry
void setRecommendedSymmetry(std::string val)
Sets the ProSHADE detected symmetry type.
Definition: ProSHADE.cpp:1231
ProSHADE_settings::setMaxSymmetryFold
void setMaxSymmetryFold(proshade_unsign maxFold)
Sets the maximum symmetry fold (well, the maximum prime symmetry fold).
Definition: ProSHADE.cpp:1403
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:83
ProSHADE_settings::setNegativeDensity
void setNegativeDensity(bool nDens)
Sets the internal variable deciding whether input files negative density should be removed.
Definition: ProSHADE.cpp:1457
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:75
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:931
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:94
ProSHADE_settings::setCorrelationMasking
void setCorrelationMasking(bool corMask)
Sets the requested map masking type in the appropriate variable.
Definition: ProSHADE.cpp:545
ProSHADE_settings::setProgressiveSphereMapping
void setProgressiveSphereMapping(bool progSphMap)
Sets the requested sphere mapping value settings approach in the appropriate variable.
Definition: ProSHADE.cpp:951
ProSHADE_settings::changeMapResolution
bool changeMapResolution
Should maps be re-sampled to obtain the required resolution?
Definition: ProSHADE_settings.hpp:51
ProSHADE_settings::setIntegrationOrder
void setIntegrationOrder(proshade_unsign intOrd)
Sets the requested order for the Gauss-Legendre integration in the appropriate variable.
Definition: ProSHADE.cpp:910
ProSHADE_settings::setRotationFunctionComputation
void setRotationFunctionComputation(bool rotfVal)
Sets whether the rotation function distance descriptor should be computed.
Definition: ProSHADE.cpp:1014
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:59
ProSHADE_settings::setOverlaySaveFile
void setOverlaySaveFile(std::string filename)
Sets the filename to which the overlay structure is to be save into.
Definition: ProSHADE.cpp:1349
ProSHADE_settings::overlayStructureName
std::string overlayStructureName
The filename to which the rotated and translated moving structure is to be saved.
Definition: ProSHADE_settings.hpp:138
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:99
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:117
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:121
ProSHADE_settings::setFourierWeightsFilename
void setFourierWeightsFilename(std::string fWgFln)
Sets the filename of the mask data that should be applied to the input map.
Definition: ProSHADE.cpp:667
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:82
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:67
ProSHADE_settings::setGroupingSmoothingFactor
void setGroupingSmoothingFactor(proshade_double smFact)
Sets the grouping smoothing factor into the proper variable.
Definition: ProSHADE.cpp:1123
ProSHADE_settings::setRequestedSymmetry
void setRequestedSymmetry(std::string val)
Sets the user requested symmetry type.
Definition: ProSHADE.cpp:1274
ProSHADE_settings::setExtraSpace
void setExtraSpace(proshade_single exSpace)
Sets the requested map extra space value in the appropriate variable.
Definition: ProSHADE.cpp:850
ProSHADE_settings::forceBounds
proshade_signed * forceBounds
These will be the boundaries to be forced upon the map.
Definition: ProSHADE_settings.hpp:96
ProSHADE_settings::fourierWeightsFileName
std::string fourierWeightsFileName
The filename from which Fourier weights data will be read from.
Definition: ProSHADE_settings.hpp:89
ProSHADE_settings::recommendedSymmetryFold
proshade_unsign recommendedSymmetryFold
The fold of the recommended symmetry C or D type, 0 otherwise.
Definition: ProSHADE_settings.hpp:129
ProSHADE_settings::inputFiles
std::vector< std::string > inputFiles
This vector contains the filenames of all input structure files.
Definition: ProSHADE_settings.hpp:43
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:810
ProSHADE_settings::setMaskSaving
void setMaskSaving(bool savMsk)
Sets whether the mask should be saved.
Definition: ProSHADE.cpp:607
ProSHADE_settings::firstModelOnly
bool firstModelOnly
Shoud only the first PDB model be used, or should all models be used?
Definition: ProSHADE_settings.hpp:46
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:728
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:125
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:95
ProSHADE_settings::setMasking
void setMasking(bool mask)
Sets the requested map masking decision value in the appropriate variable.
Definition: ProSHADE.cpp:524
ProSHADE_settings::pdbBFactorNewVal
proshade_double pdbBFactorNewVal
Change all PDB B-factors to this value (for smooth maps).
Definition: ProSHADE_settings.hpp:55
ProSHADE_settings::setRecommendedFold
void setRecommendedFold(proshade_unsign val)
Sets the ProSHADE detected symmetry fold.
Definition: ProSHADE.cpp:1254
ProSHADE_settings::removeWaters
bool removeWaters
Should all waters be removed from input PDB files?
Definition: ProSHADE_settings.hpp:45
ProSHADE_settings::setResolution
void setResolution(proshade_single resolution)
Sets the requested resolution in the appropriate variable.
Definition: ProSHADE.cpp:383
ProSHADE_settings::printSettings
void printSettings(void)
This function prints the current values in the settings object.
Definition: ProSHADE.cpp:2469
ProSHADE_run::getNoStructures
proshade_unsign getNoStructures(void)
This function returns the number of structures used.
Definition: ProSHADE.cpp:2782
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:112
ProSHADE_settings::setNormalisation
void setNormalisation(bool normalise)
Sets the requested map normalisation value in the appropriate variable.
Definition: ProSHADE.cpp:423
ProSHADE_settings::setEnergyLevelsComputation
void setEnergyLevelsComputation(bool enLevDesc)
Sets whether the energy level distance descriptor should be computed.
Definition: ProSHADE.cpp:972
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:111