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