ProSHADE  0.7.5.1 (JAN 2021)
Protein Shape Detection
pyProSHADE_data.cpp
1 
21 //==================================================== Include PyBind11 header
22 #include <pybind11/pybind11.h>
23 #include <pybind11/stl.h>
24 #include <pybind11/numpy.h>
25 
26 //==================================================== Add the ProSHADE_settings and ProSHADE_run classes to the PyBind11 module
27 void add_dataClass ( pybind11::module& pyProSHADE )
28 {
29  //================================================ Export the ProSHADE_settings class
30  pybind11::class_ < ProSHADE_internal_data::ProSHADE_data > ( pyProSHADE, "ProSHADE_data" )
31 
32  //============================================ Constructors (destructors do not need wrappers???)
33  .def ( pybind11::init < ProSHADE_settings* > ( ), pybind11::arg ( "settings" ) )
34  .def ( pybind11::init ( [] ( ProSHADE_settings* settings, std::string strName, pybind11::array_t < float, pybind11::array::c_style | pybind11::array::forcecast > mapData, proshade_single xDmSz, proshade_single yDmSz, proshade_single zDmSz, proshade_unsign xDmInd, proshade_unsign yDmInd, proshade_unsign zDmInd, proshade_signed xFr, proshade_signed yFr, proshade_signed zFr, proshade_signed xT, proshade_signed yT, proshade_signed zT, proshade_unsign inputO )
35  {
36  //== Find the array size
37  pybind11::buffer_info buf = mapData.request();
38  proshade_unsign len = buf.size;
39 
40  //== Allocate memory
41  double* npVals = new double[static_cast<proshade_unsign> ( len )];
42  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
43 
44  //== Copy from numpy to ProSHADE (I do not want to be directly changing the python memory, that screams trouble)
45  if ( buf.ndim == 1 )
46  {
47  for ( proshade_unsign iter = 0; iter < len; iter++ ) { npVals[iter] = static_cast < double > ( mapData.at(iter) ); }
48  }
49  else if ( buf.ndim == 3 )
50  {
51  float* dataPtr = reinterpret_cast < float* > ( buf.ptr );
52  for ( proshade_unsign xIt = 0; xIt < static_cast<proshade_unsign> ( buf.shape.at(0) ); xIt++ )
53  {
54  for ( proshade_unsign yIt = 0; yIt < static_cast<proshade_unsign> ( buf.shape.at(1) ); yIt++ )
55  {
56  for ( proshade_unsign zIt = 0; zIt < static_cast<proshade_unsign> ( buf.shape.at(2) ); zIt++ )
57  {
58  npVals[zIt + buf.shape.at(2) * ( yIt + buf.shape.at(1) * xIt )] = static_cast < double > ( dataPtr[zIt + buf.shape.at(2) * ( yIt + buf.shape.at(1) * xIt )] );
59  }
60  }
61  }
62  }
63  else
64  {
65  std::cerr << "!!! ProSHADE PYTHON MODULE ERROR !!! The ProSHADE_data class constructor ( ProSHADE_settings, str, numpy.ndarray, float, float, float, ... ) only supports the third argument input array in the 1D or 3D numpy.ndarray format. The supplied array has " << buf.ndim << " dims. Terminating..." << std::endl;
66  exit ( EXIT_FAILURE );
67  }
68 
69  //== Call the ProSHADE_data constructor
70  return new ProSHADE_internal_data::ProSHADE_data ( settings,
71  strName,
72  npVals,
73  static_cast<int> ( len ),
74  xDmSz,
75  yDmSz,
76  zDmSz,
77  xDmInd,
78  yDmInd,
79  zDmInd,
80  xFr,
81  yFr,
82  zFr,
83  xT,
84  yT,
85  zT,
86  inputO );
87  } ) )
88 
89  //============================================ Data I/O functions
90  .def ( "readInStructure", &ProSHADE_internal_data::ProSHADE_data::readInStructure, "This function initialises the basic ProSHADE_data variables and reads in a single structure.", pybind11::arg ( "fname" ), pybind11::arg ( "inputO" ), pybind11::arg ( "settings" ) )
91  .def ( "writeMap", &ProSHADE_internal_data::ProSHADE_data::writeMap, "Function for writing out the internal structure representation in MRC MAP format.", pybind11::arg ( "fname" ), pybind11::arg ( "title" ) = "Created by ProSHADE and written by GEMMI", pybind11::arg ( "mode" ) = 2 )
92  .def ( "writePdb", &ProSHADE_internal_data::ProSHADE_data::writePdb, "This function writes out the PDB formatted file coresponding to the structure.", pybind11::arg ( "fname" ), pybind11::arg ( "euA" ) = 0.0, pybind11::arg ( "euB" ) = 0.0, pybind11::arg ( "euG" ) = 0.0, pybind11::arg ( "trsX" ) = 0.0, pybind11::arg ( "trsY" ) = 0.0, pybind11::arg ( "trsZ" ) = 0.0, pybind11::arg ( "firstModel" ) = true )
93  .def ( "getMap",
94  [] ( ProSHADE_internal_data::ProSHADE_data &self ) -> pybind11::array_t < proshade_double >
95  {
96  //== Copy the values
97  pybind11::array_t < proshade_double > retArr = pybind11::array_t < proshade_double > ( { self.xDimIndices, self.yDimIndices, self.zDimIndices }, // Shape
98  { self.yDimIndices * self.zDimIndices * sizeof(proshade_double),
99  self.zDimIndices * sizeof(proshade_double),
100  sizeof(proshade_double) }, // C-stype strides
101  self.internalMap ); // Data
102 
103  //== Done
104  return ( retArr );
105  } )
106 
107  //============================================ Data processing functions
108  .def ( "processInternalMap", &ProSHADE_internal_data::ProSHADE_data::processInternalMap, "This function simply clusters several map manipulating functions which should be called together. These include centering, phase removal, normalisation, adding extra space, etc.", pybind11::arg ( "settings" ) )
109  .def ( "invertMirrorMap", &ProSHADE_internal_data::ProSHADE_data::invertMirrorMap, "Function for inverting the map to its mirror image.", pybind11::arg ( "settings" ) )
110  .def ( "normaliseMap", &ProSHADE_internal_data::ProSHADE_data::normaliseMap, "Function for normalising the map values to mean 0 and sd 1.", pybind11::arg ( "settings" ) )
111  .def ( "maskMap", &ProSHADE_internal_data::ProSHADE_data::maskMap, "Function for computing the map mask using blurring and X IQRs from median.", pybind11::arg ( "settings" ) )
112  .def ( "reSampleMap", &ProSHADE_internal_data::ProSHADE_data::reSampleMap, "This function changes the internal map sampling to conform to particular resolution value.", pybind11::arg ( "settings" ) )
113  .def ( "centreMapOnCOM", &ProSHADE_internal_data::ProSHADE_data::centreMapOnCOM, "This function shits the map so that its COM is in the centre of the map.", pybind11::arg ( "settings" ) )
114  .def ( "addExtraSpace", &ProSHADE_internal_data::ProSHADE_data::addExtraSpace, "This function increases the size of the map so that it can add empty space around it.", pybind11::arg ( "settings" ) )
115  .def ( "removePhaseInormation", &ProSHADE_internal_data::ProSHADE_data::removePhaseInormation, "This function removes phase from the map, effectively converting it to Patterson map.", pybind11::arg ( "settings" ) )
116  .def ( "getReBoxBoundaries",
117  [] ( ProSHADE_internal_data::ProSHADE_data &self,ProSHADE_settings* settings ) -> pybind11::array_t < proshade_signed >
118  {
119  //== Allocate output memory
120  proshade_signed* retVals = new proshade_signed[6];
121  ProSHADE_internal_misc::checkMemoryAllocation ( retVals, __FILE__, __LINE__, __func__ );
122 
123  //== If same bounds as first one are required, test if possible and return these instead
124  if ( settings->useSameBounds && ( self.inputOrder != 0 ) )
125  {
126  for ( proshade_unsign iter = 0; iter < 6; iter++ ) { retVals[iter] = settings->forceBounds[iter]; }
127  }
128  //== In this case, bounds need to be found de novo
129  else
130  {
131  //== Find the non-zero bounds
132  ProSHADE_internal_mapManip::getNonZeroBounds ( self.internalMap, self.xDimIndices, self.yDimIndices, self.zDimIndices,
133  self.xDimSize, self.yDimSize, self.zDimSize, retVals );
134 
135  //== Add the extra space
136  ProSHADE_internal_mapManip::addExtraBoundSpace ( self.xDimIndices, self.yDimIndices, self.zDimIndices,
137  self.xDimSize, self.yDimSize, self.zDimSize, retVals, settings->boundsExtraSpace );
138 
139  //== Beautify boundaries
140  ProSHADE_internal_mapManip::beautifyBoundaries ( retVals, self.xDimIndices, self.yDimIndices, self.zDimIndices, settings->boundsSimilarityThreshold, settings->verbose );
141 
142  //== If need be, save boundaries to be used for all other structure
143  if ( settings->useSameBounds && ( self.inputOrder == 0 ) ) { for ( proshade_unsign iter = 0; iter < 6; iter++ ) { settings->forceBounds[iter] = retVals[iter]; } }
144  }
145 
146  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
147  pybind11::capsule pyCapsuleReBox2 ( retVals, []( void *f ) { proshade_signed* foo = reinterpret_cast< proshade_signed* > ( f ); delete foo; } );
148 
149  //== Copy the value
150  pybind11::array_t < proshade_signed > retArr = pybind11::array_t < proshade_signed > ( { 6 }, // Shape
151  { sizeof(proshade_signed) }, // C-stype strides
152  retVals, // Data
153  pyCapsuleReBox2 ); // Capsule
154 
155  //== Done
156  return ( retArr );
157  }, "This function finds the boundaries enclosing positive map values and adds some extra space." )
158  .def ( "createNewMapFromBounds",
159  [] ( ProSHADE_internal_data::ProSHADE_data &self, pybind11::array_t < proshade_signed > bounds, ProSHADE_internal_data::ProSHADE_data* newStr, ProSHADE_settings* settings )
160  {
161  //== Allocate memory for bounds copy
162  proshade_signed* newBounds = new proshade_signed[6];
163  ProSHADE_internal_misc::checkMemoryAllocation ( newBounds, __FILE__, __LINE__, __func__ );
164 
165  //== Copy to the pointer
166  for ( proshade_unsign iter = 0; iter < 6; iter++ ) { newBounds[iter] = bounds.at(iter); }
167 
168  //== Run C++ function
169  self.createNewMapFromBounds ( settings, newStr, newBounds );
170 
171  //== Done
172  return ;
173  }, "This function creates a new structure from the calling structure and new bounds values." )
174 
175  //============================================ Data sphere mapping functions
176  .def ( "mapToSpheres", &ProSHADE_internal_data::ProSHADE_data::mapToSpheres, "This function converts the internal map onto a set of concentric spheres.", pybind11::arg ( "settings" ) )
177  .def ( "getSpherePositions", &ProSHADE_internal_data::ProSHADE_data::getSpherePositions, "This function determines the sphere positions (radii) for sphere mapping.", pybind11::arg ( "settings" ) )
178  .def ( "computeSphericalHarmonics", &ProSHADE_internal_data::ProSHADE_data::computeSphericalHarmonics, "This function computes the spherical harmonics decomposition for the whole structure.", pybind11::arg ( "settings" ) )
179 
180  //============================================ Accessor functions
181  .def ( "getXDimSize", &ProSHADE_internal_data::ProSHADE_data::getXDimSize, "This function allows access to the map size in angstroms along the X axis." )
182  .def ( "getYDimSize", &ProSHADE_internal_data::ProSHADE_data::getYDimSize, "This function allows access to the map size in angstroms along the Y axis." )
183  .def ( "getZDimSize", &ProSHADE_internal_data::ProSHADE_data::getZDimSize, "This function allows access to the map size in angstroms along the Z axis." )
184  .def ( "getXDim", &ProSHADE_internal_data::ProSHADE_data::getXDim, "This function allows access to the map size in indices along the X axis." )
185  .def ( "getYDim", &ProSHADE_internal_data::ProSHADE_data::getYDim, "This function allows access to the map size in indices along the Y axis." )
186  .def ( "getZDim", &ProSHADE_internal_data::ProSHADE_data::getZDim, "This function allows access to the map size in indices along the Z axis." )
187 
188  //============================================ Symmetry related functions
189  .def ( "computeRotationFunction", &ProSHADE_internal_data::ProSHADE_data::computeRotationFunction, "This function computes the self-rotation function for this structure and stores it internally in the ProSHADE_data object.", pybind11::arg ( "settings" ) )
190  .def ( "detectSymmetryInStructure",
192  {
193  //== Call the appropriate C++ function
194  if ( settings->usePeakSearchInRotationFunctionSpace )
195  {
196  //== Detect point groups in the angle-axis space
197  self.detectSymmetryFromAngleAxisSpace ( settings, &settings->detectedSymmetry, &settings->allDetectedCAxes );
198  }
199  else
200  {
201  //== Detect symmetry using the peak detection in rotation function space
202  self.detectSymmetryInStructure ( settings, &settings->detectedSymmetry, &settings->allDetectedCAxes );
203  }
204  }, "This function runs the symmetry detection algorithms on this structure and saves the results in the settings object.", pybind11::arg ( "settings" ) )
205  .def ( "getRecommendedSymmetryType", &ProSHADE_internal_data::ProSHADE_data::getRecommendedSymmetryType, "This function simply returns the detected recommended symmetry type.", pybind11::arg ( "settings" ) )
206  .def ( "getRecommendedSymmetryFold", &ProSHADE_internal_data::ProSHADE_data::getRecommendedSymmetryFold, "This function simply returns the detected recommended symmetry fold.", pybind11::arg ( "settings" ) )
207  .def ( "getRecommendedSymmetryAxes",
208  [] ( ProSHADE_internal_data::ProSHADE_data &self, ProSHADE_settings* settings ) -> pybind11::array_t < float >
209  {
210  //== Allocate memory for the numpy values
211  float* npVals = new float[static_cast<unsigned int> ( settings->detectedSymmetry.size() * 6 )];
212  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
213 
214  //== Copy values
215  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( settings->detectedSymmetry.size() ); iter++ ) { for ( proshade_unsign it = 0; it < 6; it++ ) { npVals[(iter*6)+it] = settings->detectedSymmetry.at(iter)[it]; } }
216 
217  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
218  pybind11::capsule pyCapsuleStrRecSym ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
219 
220  //== Copy the value
221  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<int> ( settings->detectedSymmetry.size() ), static_cast<int> ( 6 ) }, // Shape
222  { 6 * sizeof(float), sizeof(float) }, // C-stype strides
223  npVals, // Data
224  pyCapsuleStrRecSym ); // Capsule
225 
226  //== Done
227  return ( retArr );
228  }, "This function returns the recommended symmetry axes as a 2D numpy array." )
229  .def ( "getAllCSyms",
230  [] ( ProSHADE_internal_data::ProSHADE_data &self, ProSHADE_settings* settings ) -> pybind11::array_t < float >
231  {
232  //== Allocate memory for the numpy values
233  float* npVals = new float[static_cast<unsigned int> ( settings->allDetectedCAxes.size() * 6 )];
234  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
235 
236  //== Copy values
237  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( settings->allDetectedCAxes.size() ); iter++ ) { for ( proshade_unsign it = 0; it < 6; it++ ) { npVals[(iter*6)+it] = settings->allDetectedCAxes.at(iter).at(it); } }
238 
239  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
240  pybind11::capsule pyCapsuleStrSymList ( 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<int> ( settings->allDetectedCAxes.size() ), 6 }, // Shape
244  { 6 * sizeof(float), sizeof(float) }, // C-stype strides
245  npVals, // Data
246  pyCapsuleStrSymList ); // Capsule
247 
248  //== Done
249  return ( retArr );
250  }, "This function returns all symmetry axes as a 2D numpy array." )
251  .def ( "getNonCSymmetryAxesIndices",
252  [] ( ProSHADE_internal_data::ProSHADE_data &self, ProSHADE_settings* settings ) -> pybind11::dict
253  {
254  //== Initialise local variables
255  pybind11::dict retDict;
256  pybind11::list dList, tList, oList, iList;
257 
258  //== Fill in the D list
259  for ( proshade_unsign dIt = 0; dIt < static_cast<proshade_unsign> ( settings->allDetectedDAxes.size() ); dIt++ )
260  {
261  pybind11::list memList;
262  for ( proshade_unsign memIt = 0; memIt < static_cast<proshade_unsign> ( settings->allDetectedDAxes.at(dIt).size() ); memIt++ )
263  {
264  memList.append ( settings->allDetectedDAxes.at(dIt).at(memIt) );
265  }
266  dList.append ( memList );
267  }
268 
269  //== Fill in the T list
270  for ( proshade_unsign tIt = 0; tIt < static_cast<proshade_unsign> ( settings->allDetectedTAxes.size() ); tIt++ )
271  {
272  tList.append ( settings->allDetectedTAxes.at ( tIt ) );
273  }
274 
275  //== Fill in the O list
276  for ( proshade_unsign oIt = 0; oIt < static_cast<proshade_unsign> ( settings->allDetectedOAxes.size() ); oIt++ )
277  {
278  oList.append ( settings->allDetectedOAxes.at ( oIt ) );
279  }
280 
281  //== Fill in the T list
282  for ( proshade_unsign iIt = 0; iIt < static_cast<proshade_unsign> ( settings->allDetectedIAxes.size() ); iIt++ )
283  {
284  iList.append ( settings->allDetectedIAxes.at ( iIt ) );
285  }
286 
287  //== Save the lists to the dict
288  retDict[ pybind11::handle ( pybind11::str ( "D" ).ptr ( ) ) ] = dList;
289  retDict[ pybind11::handle ( pybind11::str ( "T" ).ptr ( ) ) ] = tList;
290  retDict[ pybind11::handle ( pybind11::str ( "O" ).ptr ( ) ) ] = oList;
291  retDict[ pybind11::handle ( pybind11::str ( "I" ).ptr ( ) ) ] = iList;
292 
293  //== Done
294  return ( retDict );
295  }, "This function returns array of non-C axes indices." )
296  .def ( "getAllGroupElements",
297  [] ( ProSHADE_internal_data::ProSHADE_data &self, ProSHADE_settings* settings, pybind11::array_t < int > axList, std::string groupType, proshade_double matrixTolerance ) -> pybind11::list
298  {
299  //== Sanity check
300  pybind11::buffer_info buf = axList.request();
301  if ( buf.ndim != 1 ) { std::cerr << "!!! ProSHADE PYTHON MODULE ERROR !!! The second argument to getAllGroupElements() must be a 1D numpy array stating the indices of the axes forming the group!" << std::endl; exit ( EXIT_FAILURE ); }
302 
303  //== Convert to vector of unsigns
304  std::vector< proshade_unsign > axesList;
305  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( axList.size() ); iter++ ) { ProSHADE_internal_misc::addToUnsignVector ( &axesList, axList.at(iter) ); }
306 
307  //== Get the results
308  std::vector < std::vector< proshade_double > > vals = self.getAllGroupElements ( settings, axesList, groupType, matrixTolerance );
309 
310  //== Initialise return variable
311  pybind11::list retList;
312 
313  //== Copy values
314  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ )
315  {
316  //== Allocate memory for the numpy values
317  float* npVals = new float[static_cast<unsigned int> ( 9 )];
318  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
319 
320  //== Copy values to memory
321  for ( proshade_unsign it = 0; it < 9; it++ ) { npVals[it] = vals.at(iter).at(it); }
322 
323  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
324  pybind11::capsule pyCapsuleGrpEl ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
325 
326  //== Copy the value
327  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { 3, 3 }, // Shape
328  { 3 * sizeof(float), sizeof(float) }, // C-stype strides
329  npVals, // Data
330  pyCapsuleGrpEl ); // Capsule
331 
332  //== Save matrix
333  retList.append ( retArr );
334  }
335 
336  //== Done
337  return ( retList );
338  }, "This function returns the group elements as rotation matrices of any point group described by the detected axes.", pybind11::arg ( "settings" ), pybind11::arg ( "axList" ), pybind11::arg ( "groupType" ) = "", pybind11::arg( "matrixTolerance" ) = 0.05 )
339 
340  //============================================ Overlay related functions
341  .def ( "getOverlayRotationFunction", &ProSHADE_internal_data::ProSHADE_data::getOverlayRotationFunction, "This function computes the overlay rotation function (i.e. the correlation function in SO(3) space).", pybind11::arg ( "settings" ), pybind11::arg ( "obj2" ) )
342  .def ( "getBestRotationMapPeaksEulerAngles",
343  [] ( ProSHADE_internal_data::ProSHADE_data &self, ProSHADE_settings* settings ) -> pybind11::array_t < float >
344  {
345  //== Get values
346  std::vector< proshade_double > vals = self.getBestRotationMapPeaksEulerAngles ( settings );
347 
348  //== Allocate memory for the numpy values
349  float* npVals = new float[static_cast<unsigned int> (vals.size())];
350  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
351 
352  //== Copy values
353  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] = vals.at(iter); }
354 
355  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
356  pybind11::capsule pyCapsuleEuPeak ( npVals, []( void *f ) { float* foo = reinterpret_cast< float* > ( f ); delete foo; } );
357 
358  //== Copy the value
359  pybind11::array_t < float > retArr = pybind11::array_t<float> ( { static_cast<unsigned int> (vals.size()) }, // Shape
360  { sizeof(float) }, // C-stype strides
361  npVals, // Data
362  pyCapsuleEuPeak ); // Capsule (C++ destructor, basically)
363 
364  //== Done
365  return ( retArr );
366  }, "This function returns a vector of three floats, the three Euler angles of the best peak in the rotation map.", pybind11::arg ( "settings" ) )
367  .def ( "getBestRotationMapPeaksRotationMatrix",
368  [] ( ProSHADE_internal_data::ProSHADE_data &self, ProSHADE_settings* settings ) -> pybind11::array_t < proshade_double >
369  {
370  //== Get values
371  std::vector< proshade_double > vals = self.getBestRotationMapPeaksEulerAngles ( settings );
372 
373  //== Convert Euler ZXZ to matrix
374  proshade_double* retMat = new proshade_double[9];
375  ProSHADE_internal_misc::checkMemoryAllocation ( retMat, __FILE__, __LINE__, __func__ );
376  ProSHADE_internal_maths::getRotationMatrixFromEulerZXZAngles ( vals.at(0), vals.at(1), vals.at(2), retMat );
377 
378  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
379  pybind11::capsule pyCapsuleRMPeak ( retMat, []( void *f ) { proshade_double* foo = reinterpret_cast< proshade_double* > ( f ); delete foo; } );
380 
381  //== Copy the value
382  pybind11::array_t < proshade_double > retArr = pybind11::array_t<proshade_double> ( { 3, 3 }, // Shape
383  { 3 * sizeof(proshade_double), sizeof(proshade_double) }, // C-stype strides
384  retMat, // Data
385  pyCapsuleRMPeak ); // Capsule
386 
387  //== Done
388  return ( retArr );
389  }, "This function returns a rotation matrix representing the best peak in the rotation map.", pybind11::arg ( "settings" ) )
390  .def ( "rotateMap", &ProSHADE_internal_data::ProSHADE_data::rotateMap, "This function rotates a map based on the given Euler angles.", pybind11::arg ( "settings" ), pybind11::arg ( "eulerAlpha" ), pybind11::arg ( "eulerBeta" ), pybind11::arg ( "eulerGamma" ) )
391  .def ( "zeroPaddToDims", &ProSHADE_internal_data::ProSHADE_data::zeroPaddToDims, "This function changes the size of a structure to fit the supplied new limits.", pybind11::arg ( "xDimMax" ), pybind11::arg ( "yDimMax" ), pybind11::arg ( "zDimMax" ) )
392  .def ( "computeTranslationMap", &ProSHADE_internal_data::ProSHADE_data::computeTranslationMap, "This function does the computation of the translation map and saves results internally.", pybind11::arg ( "staticStructure" ) )
393  .def ( "getOverlayTranslations",
394  [] ( ProSHADE_internal_data::ProSHADE_data &self, ProSHADE_internal_data::ProSHADE_data* staticStructure , proshade_double eulA, proshade_double eulB, proshade_double eulG) -> pybind11::dict
395  {
396  //== Get values
397  std::vector< proshade_double > vals = self.getBestTranslationMapPeaksAngstrom ( staticStructure, eulA, eulB, eulG );
398 
399  //== Initialise variables
400  pybind11::dict retDict;
401  pybind11::list rotCen, toOverlay;
402 
403  //== Copy values
404  rotCen.append ( self.originalPdbRotCenX );
405  rotCen.append ( self.originalPdbRotCenY );
406  rotCen.append ( self.originalPdbRotCenZ );
407 
408  toOverlay.append ( self.originalPdbTransX );
409  toOverlay.append ( self.originalPdbTransY );
410  toOverlay.append ( self.originalPdbTransZ );
411 
412  //== Save results to return dict
413  retDict[ pybind11::handle ( pybind11::str ( "centreOfRotation" ).ptr ( ) ) ] = rotCen;
414  retDict[ pybind11::handle ( pybind11::str ( "rotCenToOverlay" ).ptr ( ) ) ] = toOverlay;
415 
416  //== Done
417  return ( retDict );
418  }, "This function returns the vector from optimal rotation centre to origin and the optimal overlay translation vector. These two vectors allow overlaying the inputs (see documentation for details on how the two vectors should be used).", pybind11::arg ( "staticStructure" ), pybind11::arg ( "eulA" ), pybind11::arg ( "eulB" ), pybind11::arg ( "eulG" ) )
419  .def ( "translateMap", &ProSHADE_internal_data::ProSHADE_data::translateMap, "This function translates the map by a given number of Angstroms along the three axes. Please note the translation happens firstly to the whole map box and only the translation remainder that cannot be achieved by moving the box will be corrected for using reciprocal space translation within the box.", pybind11::arg ( "settings" ), pybind11::arg ( "trsX" ), pybind11::arg ( "trsY" ), pybind11::arg ( "trsZ" ) )
420 
421  //============================================ Internal arrays access functions
422  .def ( "findSHIndex",
423  [] ( ProSHADE_internal_data::ProSHADE_data &self, proshade_signed shell, proshade_signed band, proshade_signed order ) -> proshade_signed
424  {
425  //== Get value
426  proshade_signed index = seanindex ( order, band, self.spheres[shell]->getLocalBandwidth() );
427 
428  //== Done
429  return ( index );
430  }, "This function finds the correct index for given shell, band and order in the spherical harmonics array. Please note that the order is expected in range -band <= 0 <- band and NOT from 0 to ( 2 * band ) + 1." )
431  .def ( "getSphericalHarmonics",
432  [] ( ProSHADE_internal_data::ProSHADE_data &self ) -> pybind11::array_t < std::complex<proshade_double> >
433  {
434  //== Allocate memory for the numpy values
435  std::complex<proshade_double>* npVals = new std::complex<proshade_double>[static_cast<proshade_unsign> ( self.noSpheres * pow( self.maxShellBand, 2.0 ) )];
436  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
437 
438  //== Fill with zeroes as the matrix will be sparse
439  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( self.noSpheres * pow( self.maxShellBand, 2.0 ) ); iter++ ) { npVals[iter].real ( 0.0 ); npVals[iter].imag ( 0.0 ); }
440 
441  //== Initialise variables
442  proshade_signed pyPosSH;
443  proshade_signed pyPos;
444 
445  //== Copy data to new memory
446  for ( proshade_signed shIt = 0; shIt < static_cast<proshade_signed> ( self.noSpheres ); shIt++ )
447  {
448  for ( proshade_signed bnd = 0; bnd < static_cast<proshade_signed> ( self.spheres[shIt]->getLocalBandwidth() ); bnd++ )
449  {
450  for ( proshade_signed order = -bnd; order <= bnd; order++ )
451  {
452  pyPosSH = ( shIt * pow( self.maxShellBand, 2.0 ) );
453  pyPos = seanindex ( order, bnd, self.spheres[shIt]->getLocalBandwidth() );
454  npVals[pyPosSH+pyPos].real ( self.sphericalHarmonics[shIt][pyPos][0] );
455  npVals[pyPosSH+pyPos].imag ( self.sphericalHarmonics[shIt][pyPos][1] );
456  }
457  }
458  }
459 
460  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
461  pybind11::capsule pyCapsuleSHs ( npVals, []( void *f ) { std::complex<proshade_double>* foo = reinterpret_cast< std::complex<proshade_double>* > ( f ); delete foo; } );
462 
463  //== Copy the value
464  pybind11::array_t < std::complex<proshade_double> > retArr = pybind11::array_t < std::complex<proshade_double > > (
465  { static_cast<int> ( self.noSpheres ), static_cast<int> ( pow ( self.maxShellBand, 2.0 ) ) },
466  { sizeof ( std::complex < proshade_double > ) * static_cast<int> ( pow ( self.maxShellBand, 2.0 ) ), sizeof ( std::complex < proshade_double > ) },
467  npVals,
468  pyCapsuleSHs );
469 
470  //== Done
471  return ( retArr );
472  }, "This function returns a 2D numpy array of complex numbers representing the spherical harmonics computed for the structure. The first dimension of the array is the spheres (i.e. each sphere has its own array) and the second dimension is the band and order combination as given by the findSHIndex() function. Please note that while each sphere can have different number of spherical harmonics coefficients (depending on the settings.progressiveSphereMapping value), this array uses maxShellBand**2 values to make sure the length are equal for all spheres. To avoid issues, please use the findSHIndex() to correctly find the index of a particular shell, band, order combination." )
473  .def ( "getEMatrix",
474  [] ( ProSHADE_internal_data::ProSHADE_data &self ) -> pybind11::array_t < std::complex < proshade_double > >
475  {
476  //== Allocate memory for the numpy values
477  std::complex<proshade_double>* npVals = new std::complex < proshade_double >[static_cast<proshade_unsign> ( self.maxShellBand * ( ( self.maxShellBand * 2 ) + 1 ) * ( ( self.maxShellBand * 2 ) + 1 ) )];
478  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
479 
480  //== Allocate local variables
481  proshade_signed index = 0;
482 
483  //== Fill with zeroes as the matrix will be sparse
484  for ( proshade_unsign iter = 0; iter < static_cast < proshade_unsign > ( self.maxShellBand * ( ( self.maxShellBand * 2 ) + 1 ) * ( ( self.maxShellBand * 2 ) + 1 ) ); iter++ ) { npVals[iter].real ( 0.0 ); npVals[iter].imag ( 0.0 ); }
485 
486  //== Copy data to new memory
487  for ( proshade_signed bandIter = 0; bandIter < static_cast< proshade_signed > ( self.maxShellBand ); bandIter++ )
488  {
489  for ( proshade_unsign order1 = 0; order1 < ( ( bandIter * 2 ) + 1 ); order1++ )
490  {
491  for ( proshade_unsign order2 = 0; order2 < ( ( bandIter * 2 ) + 1 ); order2++ )
492  {
493  index = order2 + ( ( self.maxShellBand * 2 ) + 1 ) * ( order1 + ( ( self.maxShellBand * 2 ) + 1 ) * bandIter );
494  npVals[index].real ( self.eMatrices[bandIter][order1][order2][0] );
495  npVals[index].imag ( self.eMatrices[bandIter][order1][order2][1] );
496  }
497  }
498  }
499 
500  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
501  pybind11::capsule pyCapsuleEMs ( npVals, []( void *f ) { std::complex<proshade_double>* foo = reinterpret_cast< std::complex<proshade_double>* > ( f ); delete foo; } );
502 
503  //== Create the output object
504  pybind11::array_t < std::complex<proshade_double> > retArr = pybind11::array_t < std::complex<proshade_double > > (
505  { self.maxShellBand, ( ( self.maxShellBand * 2 ) + 1 ), ( ( self.maxShellBand * 2 ) + 1 ) },
506  { sizeof ( std::complex < proshade_double > ) * static_cast<int> ( ( ( self.maxShellBand * 2 ) + 1 ) * ( ( self.maxShellBand * 2 ) + 1 ) ),
507  sizeof ( std::complex < proshade_double > ) * static_cast<int> ( ( self.maxShellBand * 2 ) + 1 ),
508  sizeof ( std::complex < proshade_double > ) },
509  npVals,
510  pyCapsuleEMs );
511 
512  //== Done
513  return ( retArr );
514  }, "This function returns the E matrix values (these are the integral over all spheres of ( c1^(l,m) * c2^(l,m') ) values) obtained when rotation function or self-rotation function are computed. The returned array has three dimensions, first being the band, second being the order1 and third being the order2. Please note that as negative indexing does not simply work, the order indexing starts from 0 - i.e. array[1][0][0] means band 1 ; order1 = -1 and order2 = -1." )
515  .def ( "getSO3Coefficients",
516  [] ( ProSHADE_internal_data::ProSHADE_data &self ) -> pybind11::array_t < std::complex < proshade_double > >
517  {
518  //== Allocate memory for the numpy values
519  std::complex<proshade_double>* npVals = new std::complex < proshade_double >[static_cast<proshade_unsign> ( self.maxShellBand * ( ( self.maxShellBand * 2 ) + 1 ) * ( ( self.maxShellBand * 2 ) + 1 ) )];
520  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
521 
522  //== Allocate local variables
523  proshade_signed index = 0;
524 
525  //== Fill with zeroes as the matrix will be sparse
526  for ( proshade_unsign iter = 0; iter < static_cast < proshade_unsign > ( self.maxShellBand * ( ( self.maxShellBand * 2 ) + 1 ) * ( ( self.maxShellBand * 2 ) + 1 ) ); iter++ ) { npVals[iter].real ( 0.0 ); npVals[iter].imag ( 0.0 ); }
527 
528  //== Copy data to new memory
529  for ( proshade_signed bandIter = 0; bandIter < static_cast< proshade_signed > ( self.maxShellBand ); bandIter++ )
530  {
531  for ( proshade_unsign order1 = 0; order1 < ( ( bandIter * 2 ) + 1 ); order1++ )
532  {
533  for ( proshade_unsign order2 = 0; order2 < ( ( bandIter * 2 ) + 1 ); order2++ )
534  {
535  index = order2 + ( ( self.maxShellBand * 2 ) + 1 ) * ( order1 + ( ( self.maxShellBand * 2 ) + 1 ) * bandIter );
536  npVals[index].real ( self.so3Coeffs[self.so3CoeffsArrayIndex ( order1 - bandIter, order2 - bandIter, bandIter )][0] );
537  npVals[index].imag ( self.so3Coeffs[self.so3CoeffsArrayIndex ( order1 - bandIter, order2 - bandIter, bandIter )][1] );
538  }
539  }
540  }
541 
542  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
543  pybind11::capsule pyCapsuleSOCoeffs ( npVals, []( void *f ) { std::complex<proshade_double>* foo = reinterpret_cast< std::complex<proshade_double>* > ( f ); delete foo; } );
544 
545  //== Create the output object
546  pybind11::array_t < std::complex<proshade_double> > retArr = pybind11::array_t < std::complex<proshade_double > > (
547  { self.maxShellBand, ( ( self.maxShellBand * 2 ) + 1 ), ( ( self.maxShellBand * 2 ) + 1 ) },
548  { sizeof ( std::complex < proshade_double > ) * static_cast<int> ( ( ( self.maxShellBand * 2 ) + 1 ) * ( ( self.maxShellBand * 2 ) + 1 ) ),
549  sizeof ( std::complex < proshade_double > ) * static_cast<int> ( ( self.maxShellBand * 2 ) + 1 ),
550  sizeof ( std::complex < proshade_double > ) },
551  npVals,
552  pyCapsuleSOCoeffs );
553 
554  //== Done
555  return ( retArr );
556  }, "This function returns the SO(3) coefficient values (these are normalised E matrix values with sign modification) obtained when rotation function or self-rotation function are computed. The returned array has three dimensions, first being the band, second being the order1 and third being the order2. Please note that as negative indexing does not simply work, the order indexing starts from 0 - i.e. array[1][0][0] means band 1 ; order1 = -1 and order2 = -1." )
557  .def ( "getRotationFunctionMap",
558  [] ( ProSHADE_internal_data::ProSHADE_data &self ) -> pybind11::array_t < std::complex < proshade_double > >
559  {
560  //== Allocate memory for the numpy values
561  std::complex<proshade_double>* npVals = new std::complex < proshade_double >[static_cast<proshade_unsign> ( ( self.maxShellBand * 2 ) * ( self.maxShellBand * 2 ) * ( self.maxShellBand * 2 ) )];
562  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
563 
564  //== Copy the values
565  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( ( self.maxShellBand * 2 ) * ( self.maxShellBand * 2 ) * ( self.maxShellBand * 2 ) ); iter++ ) { npVals[iter].real ( self.so3CoeffsInverse[iter][0] ); npVals[iter].imag ( self.so3CoeffsInverse[iter][1] ); }
566 
567  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
568  pybind11::capsule pyCapsuleRotMap ( npVals, []( void *f ) { std::complex<proshade_double>* foo = reinterpret_cast< std::complex<proshade_double>* > ( f ); delete foo; } );
569 
570  //== Create python readable object with the correct memory access
571  pybind11::array_t < std::complex < proshade_double > > retArr = pybind11::array_t < std::complex < proshade_double > > (
572  { ( self.maxShellBand * 2 ), ( self.maxShellBand * 2 ), ( self.maxShellBand * 2 ) }, // Shape
573  { ( self.maxShellBand * 2 ) * ( self.maxShellBand * 2 ) * sizeof(std::complex < proshade_double >),
574  ( self.maxShellBand * 2 ) * sizeof(std::complex < proshade_double >),
575  sizeof(std::complex < proshade_double >) }, // C-stype strides
576  npVals, // Data
577  pyCapsuleRotMap ); // Capsule (destructor)
578 
579  //== Done
580  return ( retArr );
581  }, "This function returns the (self) rotation function as a three-dimensional map of complex numbers." )
582  .def ( "getRotationMatrixFromSOFTCoordinates",
583  [] ( ProSHADE_internal_data::ProSHADE_data &self, proshade_signed xPos, proshade_signed yPos, proshade_signed zPos ) -> pybind11::array_t < proshade_double >
584  {
585  //== Allocate memory for the numpy values
586  proshade_double* npVals = new proshade_double[9];
587  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
588 
589  //== Initialise local variables
590  proshade_double eulA, eulB, eulG;
591 
592  //== Compute the Euler angles from SOFT position
593  ProSHADE_internal_maths::getEulerZXZFromSOFTPosition ( self.maxShellBand, xPos, yPos, zPos, &eulA, &eulB, &eulG );
594 
595  //== Compute the rotation matrix
597 
598  //== Create capsule to make sure memory is released properly from the allocating language (C++ in this case)
599  pybind11::capsule pyCapsuleRMSoft ( npVals, []( void *f ) { proshade_double* foo = reinterpret_cast< proshade_double* > ( f ); delete foo; } );
600 
601  //== Create python readable object with the correct memory access
602  pybind11::array_t < proshade_double > retArr = pybind11::array_t < proshade_double > ( { 3, 3 },
603  { 3 * sizeof(proshade_double), sizeof(proshade_double) },
604  npVals,
605  pyCapsuleRMSoft );
606 
607  //== Done
608  return ( retArr );
609  }, "This function converts a given rotation function map position onto the corresponding rotation matrix.", pybind11::arg ( "xPos" ), pybind11::arg ( "yPos" ), pybind11::arg ( "zPos" ) )
610  .def ( "getTranslationFunctionMap",
611  [] ( ProSHADE_internal_data::ProSHADE_data &self ) -> pybind11::array_t < std::complex < proshade_double > >
612  {
613  //== Allocate memory for the numpy values
614  std::complex<proshade_double>* npVals = new std::complex < proshade_double >[static_cast<proshade_unsign> ( self.getXDim() * self.getYDim() * self.getZDim() )];
615  ProSHADE_internal_misc::checkMemoryAllocation ( npVals, __FILE__, __LINE__, __func__ );
616 
617  //== Copy the values
618  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( self.getXDim() * self.getYDim() * self.getZDim() ); iter++ ) { npVals[iter].real ( self.translationMap[iter][0] ); npVals[iter].imag ( self.translationMap[iter][1] ); }
619 
620  //== Create capsules to make sure memory is released properly from the allocating language (C++ in this case)
621  pybind11::capsule pyCapsuleTrsMap ( npVals, []( void *f ) { std::complex<proshade_double>* foo = reinterpret_cast< std::complex<proshade_double>* > ( f ); delete foo; } );
622 
623  //== Create python readable object with the correct memory access
624  pybind11::array_t < std::complex < proshade_double > > retArr = pybind11::array_t < std::complex < proshade_double > > (
625  { self.getXDim(), self.getYDim(), self.getZDim() }, // Shape
626  { self.getYDim() * self.getZDim() * sizeof(std::complex < proshade_double >),
627  self.getZDim() * sizeof(std::complex < proshade_double >),
628  sizeof(std::complex < proshade_double >) }, // C-stype strides
629  npVals, // Data
630  pyCapsuleTrsMap ); // Capsule (destructor)
631 
632  //== Done
633  return ( retArr );
634  }, "This function returns the translation function as a three-dimensional map of complex numbers." )
635 
636  //============================================ Member variables
637  .def_readwrite ( "fileName", &ProSHADE_internal_data::ProSHADE_data::fileName )
638  .def_readwrite ( "xDimSize", &ProSHADE_internal_data::ProSHADE_data::xDimSize )
639  .def_readwrite ( "yDimSize", &ProSHADE_internal_data::ProSHADE_data::yDimSize )
640  .def_readwrite ( "zDimSize", &ProSHADE_internal_data::ProSHADE_data::zDimSize )
641  .def_readwrite ( "aAngle", &ProSHADE_internal_data::ProSHADE_data::aAngle )
642  .def_readwrite ( "bAngle", &ProSHADE_internal_data::ProSHADE_data::bAngle )
643  .def_readwrite ( "cAngle", &ProSHADE_internal_data::ProSHADE_data::cAngle )
644  .def_readwrite ( "xDimIndices", &ProSHADE_internal_data::ProSHADE_data::xDimIndices )
645  .def_readwrite ( "yDimIndices", &ProSHADE_internal_data::ProSHADE_data::yDimIndices )
646  .def_readwrite ( "zDimIndices", &ProSHADE_internal_data::ProSHADE_data::zDimIndices )
647  .def_readwrite ( "xGridIndices", &ProSHADE_internal_data::ProSHADE_data::xGridIndices )
648  .def_readwrite ( "yGridIndices", &ProSHADE_internal_data::ProSHADE_data::yGridIndices )
649  .def_readwrite ( "zGridIndices", &ProSHADE_internal_data::ProSHADE_data::zGridIndices )
650  .def_readwrite ( "xAxisOrder", &ProSHADE_internal_data::ProSHADE_data::xAxisOrder )
651  .def_readwrite ( "yAxisOrder", &ProSHADE_internal_data::ProSHADE_data::yAxisOrder )
652  .def_readwrite ( "zAxisOrder", &ProSHADE_internal_data::ProSHADE_data::zAxisOrder )
653  .def_readwrite ( "xAxisOrigin", &ProSHADE_internal_data::ProSHADE_data::xAxisOrigin )
654  .def_readwrite ( "yAxisOrigin", &ProSHADE_internal_data::ProSHADE_data::yAxisOrigin )
655  .def_readwrite ( "zAxisOrigin", &ProSHADE_internal_data::ProSHADE_data::zAxisOrigin )
656  .def_readwrite ( "xCom", &ProSHADE_internal_data::ProSHADE_data::xCom )
657  .def_readwrite ( "yCom", &ProSHADE_internal_data::ProSHADE_data::yCom )
658  .def_readwrite ( "zCom", &ProSHADE_internal_data::ProSHADE_data::zCom )
659  .def_readwrite ( "xFrom", &ProSHADE_internal_data::ProSHADE_data::xFrom )
660  .def_readwrite ( "yFrom", &ProSHADE_internal_data::ProSHADE_data::yFrom )
661  .def_readwrite ( "zFrom", &ProSHADE_internal_data::ProSHADE_data::zFrom )
662  .def_readwrite ( "xTo", &ProSHADE_internal_data::ProSHADE_data::xTo )
663  .def_readwrite ( "yTo", &ProSHADE_internal_data::ProSHADE_data::yTo )
664  .def_readwrite ( "zTo", &ProSHADE_internal_data::ProSHADE_data::zTo )
665  .def_readwrite ( "spherePos", &ProSHADE_internal_data::ProSHADE_data::spherePos )
666  .def_readwrite ( "noSpheres", &ProSHADE_internal_data::ProSHADE_data::noSpheres )
667  .def_readwrite ( "maxShellBand", &ProSHADE_internal_data::ProSHADE_data::maxShellBand )
668  .def_readwrite ( "isEmpty", &ProSHADE_internal_data::ProSHADE_data::isEmpty )
669  .def_readwrite ( "inputOrder", &ProSHADE_internal_data::ProSHADE_data::inputOrder )
670 
671  //============================================ Description
672  .def ( "__repr__", [] ( const ProSHADE_internal_data::ProSHADE_data &a ) { return "<ProSHADE_data class object> (This class contains all information, results and available functionalities for a structure)"; } );
673 }
ProSHADE_internal_mapManip::addExtraBoundSpace
void addExtraBoundSpace(proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_signed *&bounds, proshade_single extraSpace)
This function takes a set of bounds and adds a given number of Angstroms to them.
Definition: ProSHADE_mapManip.cpp:1122
ProSHADE_internal_data::ProSHADE_data::computeSphericalHarmonics
void computeSphericalHarmonics(ProSHADE_settings *settings)
This function computes the spherical harmonics decomposition for the whole structure.
Definition: ProSHADE_data.cpp:1542
ProSHADE_internal_data::ProSHADE_data::getOverlayRotationFunction
void getOverlayRotationFunction(ProSHADE_settings *settings, ProSHADE_internal_data::ProSHADE_data *obj2)
This function computes the overlay rotation function (i.e. the correlation function in SO(3) space).
Definition: ProSHADE_overlay.cpp:35
ProSHADE_internal_data::ProSHADE_data::fileName
std::string fileName
This is the original file from which the data were obtained.
Definition: ProSHADE_data.hpp:52
ProSHADE_internal_data::ProSHADE_data::zFrom
proshade_signed zFrom
This is the starting index along the z axis.
Definition: ProSHADE_data.hpp:112
ProSHADE_internal_data::ProSHADE_data::xDimIndices
proshade_unsign xDimIndices
This is the size of the map cell x dimension in indices.
Definition: ProSHADE_data.hpp:65
ProSHADE_internal_data::ProSHADE_data::getZDimSize
proshade_single getZDimSize(void)
This function allows access to the map size in angstroms along the Z axis.
Definition: ProSHADE_data.cpp:3285
ProSHADE_internal_data::ProSHADE_data::xGridIndices
proshade_unsign xGridIndices
As far as I know, this is identical to the xDimIndices.
Definition: ProSHADE_data.hpp:68
ProSHADE_internal_data::ProSHADE_data::computeRotationFunction
void computeRotationFunction(ProSHADE_settings *settings)
This function computes the self-rotation function for this structure.
Definition: ProSHADE_symmetry.cpp:34
ProSHADE_internal_data::ProSHADE_data::xAxisOrder
proshade_unsign xAxisOrder
This is the order of the x axis.
Definition: ProSHADE_data.hpp:71
ProSHADE_internal_data::ProSHADE_data::getRecommendedSymmetryFold
proshade_unsign getRecommendedSymmetryFold(ProSHADE_settings *settings)
This function simply returns the detected recommended symmetry fold.
Definition: ProSHADE_data.cpp:3743
ProSHADE_internal_data::ProSHADE_data::cAngle
proshade_single cAngle
This is the angle c of the map cell in degrees.
Definition: ProSHADE_data.hpp:64
ProSHADE_internal_data::ProSHADE_data::yCom
proshade_double yCom
The COM of the map after processing along the Y-axis.
Definition: ProSHADE_data.hpp:78
ProSHADE_settings::boundsExtraSpace
proshade_single boundsExtraSpace
The number of extra angstroms to be added to all re-boxing bounds just for safety.
Definition: ProSHADE_settings.hpp:137
ProSHADE_settings::allDetectedDAxes
std::vector< std::vector< proshade_unsign > > allDetectedDAxes
The vector of all detected dihedral symmetry axes indices in allDetectedCAxes.
Definition: ProSHADE_settings.hpp:191
ProSHADE_internal_data::ProSHADE_data::getXDimSize
proshade_single getXDimSize(void)
This function allows access to the map size in angstroms along the X axis.
Definition: ProSHADE_data.cpp:3265
ProSHADE_internal_data::ProSHADE_data
This class contains all inputed and derived data for a single structure.
Definition: ProSHADE_data.hpp:49
ProSHADE_internal_data::ProSHADE_data::noSpheres
proshade_unsign noSpheres
The number of spheres with map projected onto them.
Definition: ProSHADE_data.hpp:119
ProSHADE_internal_data::ProSHADE_data::getYDimSize
proshade_single getYDimSize(void)
This function allows access to the map size in angstroms along the Y axis.
Definition: ProSHADE_data.cpp:3275
ProSHADE_internal_data::ProSHADE_data::inputOrder
proshade_unsign inputOrder
This value is the input order - it is useful to know for writing out files, so that they would not ov...
Definition: ProSHADE_data.hpp:140
ProSHADE_internal_data::ProSHADE_data::zDimSize
proshade_single zDimSize
This is the size of the map cell z dimension in Angstroms.
Definition: ProSHADE_data.hpp:61
ProSHADE_internal_data::ProSHADE_data::maxShellBand
proshade_unsign maxShellBand
The maximum band for any shell of the object.
Definition: ProSHADE_data.hpp:123
ProSHADE_internal_data::ProSHADE_data::getSpherePositions
void getSpherePositions(ProSHADE_settings *settings)
This function determines the sphere positions (radii) for sphere mapping.
Definition: ProSHADE_data.cpp:1442
ProSHADE_internal_data::ProSHADE_data::invertMirrorMap
void invertMirrorMap(ProSHADE_settings *settings)
Function for inverting the map to its mirror image.
Definition: ProSHADE_data.cpp:894
ProSHADE_internal_maths::getRotationMatrixFromEulerZXZAngles
void getRotationMatrixFromEulerZXZAngles(proshade_double eulerAlpha, proshade_double eulerBeta, proshade_double eulerGamma, proshade_double *matrix)
Function to find the rotation matrix from Euler angles (ZXZ convention).
Definition: ProSHADE_maths.cpp:1005
ProSHADE_settings::allDetectedTAxes
std::vector< proshade_unsign > allDetectedTAxes
The vector of all detected tetrahedral symmetry axes indices in allDetectedCAxes.
Definition: ProSHADE_settings.hpp:192
ProSHADE_settings::allDetectedOAxes
std::vector< proshade_unsign > allDetectedOAxes
The vector of all detected octahedral symmetry axes indices in allDetectedCAxes.
Definition: ProSHADE_settings.hpp:193
ProSHADE_internal_data::ProSHADE_data::centreMapOnCOM
void centreMapOnCOM(ProSHADE_settings *settings)
This function shits the map so that its COM is in the centre of the map.
Definition: ProSHADE_data.cpp:1228
ProSHADE_settings::verbose
proshade_signed verbose
Should the software report on the progress, or just be quiet? Value between -1 (nothing) and 4 (loud)
Definition: ProSHADE_settings.hpp:185
ProSHADE_internal_mapManip::beautifyBoundaries
void beautifyBoundaries(proshade_signed *&bounds, proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim, proshade_signed boundsDiffThres, proshade_signed verbose)
Function for modifying boundaries to a mathematically more pleasant values.
Definition: ProSHADE_mapManip.cpp:1897
ProSHADE_internal_data::ProSHADE_data::yDimIndices
proshade_unsign yDimIndices
This is the size of the map cell y dimension in indices.
Definition: ProSHADE_data.hpp:66
ProSHADE_internal_data::ProSHADE_data::yGridIndices
proshade_unsign yGridIndices
As far as I know, this is identical to the yDimIndices.
Definition: ProSHADE_data.hpp:69
ProSHADE_internal_data::ProSHADE_data::zAxisOrigin
proshade_signed zAxisOrigin
This is the origin position along the z axis.
Definition: ProSHADE_data.hpp:76
ProSHADE_internal_data::ProSHADE_data::getXDim
proshade_unsign getXDim(void)
This function allows access to the map size in indices along the X axis.
Definition: ProSHADE_data.cpp:3295
ProSHADE_internal_data::ProSHADE_data::zTo
proshade_signed zTo
This is the final index along the z axis.
Definition: ProSHADE_data.hpp:115
ProSHADE_internal_data::ProSHADE_data::xTo
proshade_signed xTo
This is the final index along the x axis.
Definition: ProSHADE_data.hpp:113
ProSHADE_internal_data::ProSHADE_data::rotateMap
void rotateMap(ProSHADE_settings *settings, proshade_double eulerAlpha, proshade_double eulerBeta, proshade_double eulerGamma)
This function rotates a map based on the given Euler angles.
Definition: ProSHADE_overlay.cpp:736
ProSHADE_internal_data::ProSHADE_data::mapToSpheres
void mapToSpheres(ProSHADE_settings *settings)
This function converts the internal map onto a set of concentric spheres.
Definition: ProSHADE_data.cpp:1499
ProSHADE_internal_data::ProSHADE_data::xAxisOrigin
proshade_signed xAxisOrigin
This is the origin position along the x axis.
Definition: ProSHADE_data.hpp:74
ProSHADE_internal_data::ProSHADE_data::yTo
proshade_signed yTo
This is the final index along the y axis.
Definition: ProSHADE_data.hpp:114
ProSHADE_internal_data::ProSHADE_data::xDimSize
proshade_single xDimSize
This is the size of the map cell x dimension in Angstroms.
Definition: ProSHADE_data.hpp:59
ProSHADE_internal_data::ProSHADE_data::yAxisOrigin
proshade_signed yAxisOrigin
This is the origin position along the y axis.
Definition: ProSHADE_data.hpp:75
ProSHADE_internal_data::ProSHADE_data::computeTranslationMap
void computeTranslationMap(ProSHADE_internal_data::ProSHADE_data *obj1)
This function does the computation of the translation map and saves results internally.
Definition: ProSHADE_overlay.cpp:389
ProSHADE_internal_data::ProSHADE_data::yAxisOrder
proshade_unsign yAxisOrder
This is the order of the y axis.
Definition: ProSHADE_data.hpp:72
ProSHADE_internal_data::ProSHADE_data::xFrom
proshade_signed xFrom
This is the starting index along the x axis.
Definition: ProSHADE_data.hpp:110
ProSHADE_internal_data::ProSHADE_data::writeMap
void writeMap(std::string fName, std::string title="Created by ProSHADE and written by GEMMI", int mode=2)
Function for writing out the internal structure representation in MRC MAP format.
Definition: ProSHADE_data.cpp:741
ProSHADE_internal_mapManip::getNonZeroBounds
void getNonZeroBounds(proshade_double *map, proshade_signed xDim, proshade_signed yDim, proshade_signed zDim, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_signed *&ret)
Function for finding the map boundaries enclosing positive only values.
Definition: ProSHADE_mapManip.cpp:1065
ProSHADE_internal_data::ProSHADE_data::getYDim
proshade_unsign getYDim(void)
This function allows access to the map size in indices along the Y axis.
Definition: ProSHADE_data.cpp:3305
ProSHADE_internal_data::ProSHADE_data::maskMap
void maskMap(ProSHADE_settings *settings)
Function for computing the map mask using blurring and X IQRs from median.
Definition: ProSHADE_data.cpp:995
ProSHADE_settings
This class stores all the settings and is passed to the executive classes instead of a multitude of p...
Definition: ProSHADE_settings.hpp:86
ProSHADE_internal_data::ProSHADE_data::xCom
proshade_double xCom
The COM of the map after processing along the X-axis.
Definition: ProSHADE_data.hpp:77
ProSHADE_internal_data::ProSHADE_data::processInternalMap
void processInternalMap(ProSHADE_settings *settings)
This function simply clusters several other functions which should be called together.
Definition: ProSHADE_data.cpp:1398
ProSHADE_internal_data::ProSHADE_data::zCom
proshade_double zCom
The COM of the map after processing along the Z-axis.
Definition: ProSHADE_data.hpp:79
ProSHADE_internal_data::ProSHADE_data::spherePos
std::vector< proshade_single > spherePos
Vector of sphere radii from the centre of the map.
Definition: ProSHADE_data.hpp:118
ProSHADE_internal_data::ProSHADE_data::getZDim
proshade_unsign getZDim(void)
This function allows access to the map size in indices along the Z axis.
Definition: ProSHADE_data.cpp:3315
ProSHADE_settings::boundsSimilarityThreshold
proshade_signed boundsSimilarityThreshold
Number of indices which can be added just to make sure same size in indices is achieved.
Definition: ProSHADE_settings.hpp:138
ProSHADE_internal_data::ProSHADE_data::writePdb
void writePdb(std::string fName, proshade_double euA=0.0, proshade_double euB=0.0, proshade_double euG=0.0, proshade_double trsX=0.0, proshade_double trsY=0.0, proshade_double trsZ=0.0, bool firstModel=true)
This function writes out the PDB formatted file coresponding to the structure so that its COM is at s...
Definition: ProSHADE_data.cpp:806
ProSHADE_internal_data::ProSHADE_data::zAxisOrder
proshade_unsign zAxisOrder
This is the order of the z axis.
Definition: ProSHADE_data.hpp:73
ProSHADE_internal_data::ProSHADE_data::reSampleMap
void reSampleMap(ProSHADE_settings *settings)
This function changes the internal map sampling to conform to particular resolution value.
Definition: ProSHADE_data.cpp:1156
ProSHADE_internal_data::ProSHADE_data::translateMap
void translateMap(ProSHADE_settings *settings, proshade_double trsX, proshade_double trsY, proshade_double trsZ)
This function simply translates the map by a given number of Angstroms along the three axes.
Definition: ProSHADE_overlay.cpp:795
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_internal_data::ProSHADE_data::bAngle
proshade_single bAngle
This is the angle b of the map cell in degrees.
Definition: ProSHADE_data.hpp:63
ProSHADE_internal_maths::getEulerZXZFromSOFTPosition
void getEulerZXZFromSOFTPosition(proshade_signed band, proshade_signed x, proshade_signed y, proshade_signed z, proshade_double *eulerAlpha, proshade_double *eulerBeta, proshade_double *eulerGamma)
Function to find Euler angles (ZXZ convention) from index position in the inverse SOFT map.
Definition: ProSHADE_maths.cpp:961
ProSHADE_internal_data::ProSHADE_data::aAngle
proshade_single aAngle
This is the angle a of the map cell in degrees.
Definition: ProSHADE_data.hpp:62
ProSHADE_internal_data::ProSHADE_data::isEmpty
bool isEmpty
This variable stated whether the class contains any information.
Definition: ProSHADE_data.hpp:139
ProSHADE_settings::allDetectedIAxes
std::vector< proshade_unsign > allDetectedIAxes
The vector of all detected icosahedral symmetry axes indices in allDetectedCAxes.
Definition: ProSHADE_settings.hpp:194
ProSHADE_settings::forceBounds
proshade_signed * forceBounds
These will be the boundaries to be forced upon the map.
Definition: ProSHADE_settings.hpp:140
ProSHADE_settings::detectedSymmetry
std::vector< proshade_double * > detectedSymmetry
The vector of detected symmetry axes.
Definition: ProSHADE_settings.hpp:189
ProSHADE_internal_data::ProSHADE_data::yFrom
proshade_signed yFrom
This is the starting index along the y axis.
Definition: ProSHADE_data.hpp:111
ProSHADE_internal_misc::addToUnsignVector
void addToUnsignVector(std::vector< proshade_unsign > *vecToAddTo, proshade_unsign elementToAdd)
Adds the element to the vector.
Definition: ProSHADE_misc.cpp:99
ProSHADE_settings::useSameBounds
bool useSameBounds
Switch to say that the same boundaries as used for the first should be used for all input maps.
Definition: ProSHADE_settings.hpp:139
ProSHADE_internal_data::ProSHADE_data::zDimIndices
proshade_unsign zDimIndices
This is the size of the map cell z dimension in indices.
Definition: ProSHADE_data.hpp:67
ProSHADE_internal_data::ProSHADE_data::getRecommendedSymmetryType
std::string getRecommendedSymmetryType(ProSHADE_settings *settings)
This function simply returns the detected recommended symmetry type.
Definition: ProSHADE_data.cpp:3732
ProSHADE_internal_data::ProSHADE_data::zGridIndices
proshade_unsign zGridIndices
As far as I know, this is identical to the zDimIndices.
Definition: ProSHADE_data.hpp:70
ProSHADE_settings::usePeakSearchInRotationFunctionSpace
bool usePeakSearchInRotationFunctionSpace
This variable switch decides whether symmetry detection will be done using peak search in rotation fu...
Definition: ProSHADE_settings.hpp:176
ProSHADE_internal_data::ProSHADE_data::normaliseMap
void normaliseMap(ProSHADE_settings *settings)
Function for normalising the map values to mean 0 and sd 1..
Definition: ProSHADE_data.cpp:948
ProSHADE_internal_data::ProSHADE_data::yDimSize
proshade_single yDimSize
This is the size of the map cell y dimension in Angstroms.
Definition: ProSHADE_data.hpp:60
ProSHADE_settings::allDetectedCAxes
std::vector< std::vector< proshade_double > > allDetectedCAxes
The vector of all detected cyclic symmetry axes.
Definition: ProSHADE_settings.hpp:190
ProSHADE_internal_data::ProSHADE_data::removePhaseInormation
void removePhaseInormation(ProSHADE_settings *settings)
This function removes phase from the map, effectively converting it to Patterson map.
Definition: ProSHADE_data.cpp:3032
ProSHADE_internal_data::ProSHADE_data::zeroPaddToDims
void zeroPaddToDims(proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim)
This function changes the size of a structure to fit the supplied new limits.
Definition: ProSHADE_overlay.cpp:598
ProSHADE_internal_data::ProSHADE_data::readInStructure
void readInStructure(std::string fName, proshade_unsign inputO, ProSHADE_settings *settings)
This function initialises the basic ProSHADE_data variables and reads in a single structure.
Definition: ProSHADE_data.cpp:447
ProSHADE_internal_data::ProSHADE_data::addExtraSpace
void addExtraSpace(ProSHADE_settings *settings)
This function increases the size of the map so that it can add empty space around it.
Definition: ProSHADE_data.cpp:1289