23 #include <pybind11/pybind11.h>
24 #include <pybind11/stl.h>
25 #include <pybind11/numpy.h>
28 void add_dataClass ( pybind11::module& pyProSHADE )
31 pybind11::class_ < ProSHADE_internal_data::ProSHADE_data > ( pyProSHADE,
"ProSHADE_data" )
34 .def ( pybind11::init ( ) )
35 .def ( pybind11::init ( [] ( 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 )
38 pybind11::buffer_info buf = mapData.request();
39 proshade_unsign len =
static_cast< proshade_unsign
> ( buf.size );
42 double* npVals =
new double[
static_cast<proshade_unsign
> ( len )];
48 for ( proshade_unsign iter = 0; iter < len; iter++ ) { npVals[iter] =
static_cast < double > ( mapData.at(iter) ); }
50 else if ( buf.ndim == 3 )
52 float* dataPtr =
reinterpret_cast < float*
> ( buf.ptr );
53 for ( proshade_unsign xIt = 0; xIt < static_cast<proshade_unsign> ( buf.shape.at(0) ); xIt++ )
55 for ( proshade_unsign yIt = 0; yIt < static_cast<proshade_unsign> ( buf.shape.at(1) ); yIt++ )
57 for ( proshade_unsign zIt = 0; zIt < static_cast<proshade_unsign> ( buf.shape.at(2) ); zIt++ )
59 npVals[zIt +
static_cast< proshade_unsign
> ( buf.shape.at(2) ) * ( yIt +
static_cast< proshade_unsign
> ( buf.shape.at(1) ) * xIt )] =
static_cast < double > ( dataPtr[zIt +
static_cast< proshade_unsign
> ( buf.shape.at(2) ) * ( yIt +
static_cast< proshade_unsign
> ( buf.shape.at(1) ) * xIt )] );
66 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;
67 exit ( EXIT_FAILURE );
73 static_cast<int> ( len ),
90 .def (
"readInStructure",
94 pybind11::buffer_info maskArr_buf = maskArr.request();
95 pybind11::buffer_info weightsArr_buf = weightsArr.request();
98 if ( ( maskArr_buf.shape.at(0) != 0 ) && ( maskArr_buf.shape.at(1) != 0 ) && ( maskArr_buf.shape.at(2) != 0 ) )
101 if ( maskArr_buf.ndim != 3 ) { std::cerr <<
"!!! ProSHADE PYTHON MODULE ERROR !!! The fourth argument to readInStructure() must be a 3D numpy array or empty and the dimensions must match!" << std::endl; exit ( EXIT_FAILURE ); }
104 proshade_double* mskArr =
new proshade_double[maskArr_buf.shape.at(0)*maskArr_buf.shape.at(1)*maskArr_buf.shape.at(2)];
108 proshade_double* arrStart =
static_cast< proshade_double*
> ( maskArr_buf.ptr );
109 for (
size_t iter = 0; iter < static_cast< size_t > ( maskArr_buf.shape.at(0)*maskArr_buf.shape.at(1)*maskArr_buf.shape.at(2) ); iter++ ) { mskArr[iter] = arrStart[iter]; }
112 if ( ( weightsArr_buf.shape.at(0) != 0 ) && ( weightsArr_buf.shape.at(1) != 0 ) && ( weightsArr_buf.shape.at(2) != 0 ) )
115 if ( weightsArr_buf.ndim != 3 ) { std::cerr <<
"!!! ProSHADE PYTHON MODULE ERROR !!! The fifth argument to readInStructure() must be a 3D numpy array or empty and the dimensions must match!" << std::endl; exit ( EXIT_FAILURE ); }
118 proshade_double* wghArr =
new proshade_double[weightsArr_buf.shape.at(0)*weightsArr_buf.shape.at(1)*weightsArr_buf.shape.at(2)];
122 proshade_double* arr2Start =
static_cast< proshade_double*
> ( weightsArr_buf.ptr );
123 for (
size_t iter = 0; iter < static_cast< size_t > ( weightsArr_buf.shape.at(0)*weightsArr_buf.shape.at(1)*weightsArr_buf.shape.at(2) ); iter++ ) { wghArr[iter] = arr2Start[iter]; }
126 self.readInStructure ( fName, inputO, settings, mskArr,
static_cast< proshade_unsign
> ( maskArr_buf.shape.at(0) ),
static_cast< proshade_unsign
> ( maskArr_buf.shape.at(1) ),
static_cast< proshade_unsign
> ( maskArr_buf.shape.at(2) ), wghArr,
static_cast< proshade_unsign
> ( weightsArr_buf.shape.at(0) ),
static_cast< proshade_unsign
> ( weightsArr_buf.shape.at(1) ),
static_cast< proshade_unsign
> ( weightsArr_buf.shape.at(2) ) );
134 self.readInStructure ( fName, inputO, settings, mskArr,
static_cast< proshade_unsign
> ( maskArr_buf.shape.at(0) ),
static_cast< proshade_unsign
> ( maskArr_buf.shape.at(1) ),
static_cast< proshade_unsign
> ( maskArr_buf.shape.at(2) ),
nullptr, 0, 0, 0 );
143 if ( ( weightsArr_buf.shape.at(0) != 0 ) && ( weightsArr_buf.shape.at(1) != 0 ) && ( weightsArr_buf.shape.at(2) != 0 ) )
146 if ( weightsArr_buf.ndim != 3 ) { std::cerr <<
"!!! ProSHADE PYTHON MODULE ERROR !!! The fifth argument to readInStructure() must be a 3D numpy array or empty and the dimensions must match!" << std::endl; exit ( EXIT_FAILURE ); }
149 proshade_double* wghArr =
new proshade_double[weightsArr_buf.shape.at(0)*weightsArr_buf.shape.at(1)*weightsArr_buf.shape.at(2)];
153 proshade_double* arr2Start =
static_cast< proshade_double*
> ( weightsArr_buf.ptr );
154 for (
size_t iter = 0; iter < static_cast< size_t > ( weightsArr_buf.shape.at(0)*weightsArr_buf.shape.at(1)*weightsArr_buf.shape.at(2) ); iter++ ) { wghArr[iter] = arr2Start[iter]; }
157 self.readInStructure ( fName, inputO, settings,
nullptr, 0, 0, 0, wghArr,
static_cast< proshade_unsign
> ( weightsArr_buf.shape.at(0) ),
static_cast< proshade_unsign
> ( weightsArr_buf.shape.at(1) ),
static_cast< proshade_unsign
> ( weightsArr_buf.shape.at(2) ) );
165 self.readInStructure ( fName, inputO, settings,
nullptr, 0, 0, 0,
nullptr, 0, 0, 0 );
171 },
"This function returns the group elements as rotation matrices of any point group described by the detected axes.", pybind11::arg (
"fName" ), pybind11::arg (
"inputO" ), pybind11::arg (
"settings" ), pybind11::arg(
"maskArr" ) = pybind11::array_t < proshade_double > (), pybind11::arg(
"weightsArr" ) = pybind11::array_t < proshade_double > () )
172 .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 )
173 .def (
"writePdb", &
ProSHADE_internal_data::ProSHADE_data::writePdb,
"This function writes out the co-ordinates file with ProSHADE type rotation and translation applied.", 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 (
"rotX" ) = 0.0, pybind11::arg (
"rotY" ) = 0.0, pybind11::arg (
"rotZ" ) = 0.0, pybind11::arg (
"firstModel" ) = true )
174 .def (
"writeGemmi", &
ProSHADE_internal_data::ProSHADE_data::writeGemmi,
"This function writes out the gemmi::Structure object with ProSHADE type rotation and translation applied.", pybind11::arg (
"fname" ), pybind11::arg (
"gemmiStruct" ), 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 (
"rotX" ) = 0.0, pybind11::arg (
"rotY" ) = 0.0, pybind11::arg (
"rotZ" ) = 0.0, pybind11::arg (
"firstModel" ) = true )
179 pybind11::array_t < proshade_double > retArr = pybind11::array_t < proshade_double > ( {
self.xDimIndices,
self.yDimIndices,
self.zDimIndices },
180 {
self.yDimIndices *
self.zDimIndices *
sizeof(proshade_double),
181 self.zDimIndices *
sizeof(proshade_double),
182 sizeof(proshade_double) },
190 .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" ) )
200 .def (
"getReBoxBoundaries",
204 proshade_signed* retVals =
new proshade_signed[6];
210 for ( proshade_unsign iter = 0; iter < 6; iter++ ) { retVals[iter] = settings->
forceBounds[iter]; }
217 static_cast< proshade_signed
> (
self.xDimIndices ),
218 static_cast< proshade_signed
> (
self.yDimIndices ),
219 static_cast< proshade_signed
> (
self.zDimIndices ),
224 self.xDimSize,
self.yDimSize,
self.zDimSize, retVals, settings->
boundsExtraSpace );
230 if ( settings->
useSameBounds && (
self.inputOrder == 0 ) ) {
for ( proshade_unsign iter = 0; iter < 6; iter++ ) { settings->
forceBounds[iter] = retVals[iter]; } }
234 pybind11::capsule pyCapsuleReBox2 ( retVals, [](
void *f ) { proshade_signed* foo =
reinterpret_cast< proshade_signed*
> ( f );
delete foo; } );
237 pybind11::array_t < proshade_signed > retArr = pybind11::array_t < proshade_signed > ( { 6 },
238 {
sizeof(proshade_signed) },
244 },
"This function finds the boundaries enclosing positive map values and adds some extra space." )
245 .def (
"createNewMapFromBounds",
249 proshade_signed* newBounds =
new proshade_signed[6];
253 for ( proshade_unsign iter = 0; iter < 6; iter++ ) { newBounds[iter] = bounds.at(iter); }
256 self.createNewMapFromBounds ( settings, newStr, newBounds );
260 },
"This function creates a new structure from the calling structure and new bounds values." )
277 .def (
"detectSymmetryInStructure",
281 self.detectSymmetryFromAngleAxisSpace ( settings );
282 },
"This function runs the symmetry detection algorithms on this structure and saves the results in the settings object.", pybind11::arg (
"settings" ) )
283 .def (
"reRunSymmetryDetectionThreshold",
287 fftw_complex* fCoeffsCut;
288 proshade_double **bindata, *fscByBin;
289 proshade_signed *cutIndices, *binCounts, noBins, cutXDim, cutYDim, cutZDim;
290 self.prepareFSCFourierMemory ( cutIndices, fCoeffsCut, &noBins, bindata, binCounts, fscByBin, settings->
requestedResolution, &cutXDim, &cutYDim, &cutZDim );
293 self.determineRecommendedSymmetry ( settings, threshold, cutIndices, fCoeffsCut, noBins, bindata, binCounts, fscByBin, cutXDim, cutYDim, cutZDim );
296 for (
size_t binIt = 0; binIt < static_cast< size_t > ( noBins ); binIt++ ) {
delete[] bindata[binIt]; }
300 fftw_free ( fCoeffsCut );
301 },
"This function runs the symmetry detection algorithms on this structure and saves the results in the settings object.", pybind11::arg (
"settings" ), pybind11::arg (
"threshold" ) )
304 .def (
"getRecommendedSymmetryAxes",
308 float* npVals =
new float[
static_cast<unsigned int> (
self.recommendedSymmetryValues.size() * 7 )];
312 for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> (
self.recommendedSymmetryValues.size() ); iter++ ) {
for ( proshade_unsign it = 0; it < 7; it++ ) { npVals[(iter*7)+it] =
static_cast< float > (
self.recommendedSymmetryValues.at(iter)[it] ); } }
315 pybind11::capsule pyCapsuleStrRecSym ( npVals, [](
void *f ) {
float* foo =
reinterpret_cast< float*
> ( f );
delete foo; } );
318 pybind11::array_t < float > retArr = pybind11::array_t<float> ( {
static_cast<int> (
self.recommendedSymmetryValues.size() ),
static_cast<int> ( 7 ) },
319 { 7 *
sizeof(float),
sizeof(
float) },
321 pyCapsuleStrRecSym );
325 },
"This function returns the recommended symmetry axes as a 2D numpy array." )
326 .def (
"getAllCSyms",
330 float* npVals =
new float[
static_cast<unsigned int> (
self.cyclicSymmetries.size() * 7 )];
334 for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> (
self.cyclicSymmetries.size() ); iter++ ) {
for ( proshade_unsign it = 0; it < 7; it++ ) { npVals[(iter*7)+it] =
static_cast< float > (
self.cyclicSymmetries.at(iter)[it] ); } }
337 pybind11::capsule pyCapsuleStrSymList ( npVals, [](
void *f ) {
float* foo =
reinterpret_cast< float*
> ( f );
delete foo; } );
340 pybind11::array_t < float > retArr = pybind11::array_t<float> ( {
static_cast<int> (
self.cyclicSymmetries.size() ), 7 },
341 { 7 *
sizeof(float),
sizeof(
float) },
343 pyCapsuleStrSymList );
347 },
"This function returns all symmetry axes as a 2D numpy array." )
348 .def (
"getNonCSymmetryAxesIndices",
352 pybind11::dict retDict;
353 pybind11::list dList, tList, oList, iList;
356 for ( proshade_unsign dIt = 0; dIt < static_cast<proshade_unsign> (
self.dihedralSymmetries.size() ); dIt++ )
358 pybind11::list memList;
359 for ( proshade_unsign memIt = 0; memIt < static_cast<proshade_unsign> (
self.dihedralSymmetries.at(dIt).size() ); memIt++ )
361 pybind11::list axList;
362 for (
size_t vIt = 0; vIt < 7; vIt++ ) { axList.append (
self.dihedralSymmetries.at(dIt).at(memIt)[vIt] ); }
363 memList.append ( axList );
365 dList.append ( memList );
369 for ( proshade_unsign tIt = 0; tIt < static_cast<proshade_unsign> (
self.tetrahedralSymmetries.size() ); tIt++ )
371 pybind11::list axList;
372 for (
size_t vIt = 0; vIt < 7; vIt++ ) { axList.append (
self.tetrahedralSymmetries.at(tIt)[vIt] ); }
373 tList.append ( axList );
377 for ( proshade_unsign oIt = 0; oIt < static_cast<proshade_unsign> (
self.octahedralSymmetries.size() ); oIt++ )
379 pybind11::list axList;
380 for (
size_t vIt = 0; vIt < 7; vIt++ ) { axList.append (
self.octahedralSymmetries.at(oIt)[vIt] ); }
381 oList.append ( axList );
385 for ( proshade_unsign iIt = 0; iIt < static_cast<proshade_unsign> (
self.icosahedralSymmetries.size() ); iIt++ )
387 pybind11::list axList;
388 for (
size_t vIt = 0; vIt < 7; vIt++ ) { axList.append (
self.icosahedralSymmetries.at(iIt)[vIt] ); }
389 iList.append ( axList );
393 retDict[ pybind11::handle ( pybind11::str (
"D" ).ptr ( ) ) ] = dList;
394 retDict[ pybind11::handle ( pybind11::str (
"T" ).ptr ( ) ) ] = tList;
395 retDict[ pybind11::handle ( pybind11::str (
"O" ).ptr ( ) ) ] = oList;
396 retDict[ pybind11::handle ( pybind11::str (
"I" ).ptr ( ) ) ] = iList;
400 },
"This function returns array of non-C axes indices." )
401 .def (
"getAllGroupElements",
405 pybind11::buffer_info buf = axList.request();
406 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 ); }
409 std::vector< proshade_unsign > axesList;
410 for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( axList.size() ); iter++ ) {
ProSHADE_internal_misc::addToUnsignVector ( &axesList,
static_cast< proshade_unsign
> ( axList.at(iter) ) ); }
413 std::vector < std::vector< proshade_double > > vals =
self.getAllGroupElements ( axesList, groupType, matrixTolerance );
416 pybind11::list retList;
419 for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ )
422 float* npVals =
new float[
static_cast<unsigned int> ( 9 )];
426 for ( proshade_unsign it = 0; it < 9; it++ ) { npVals[it] =
static_cast< float > ( vals.at(iter).at(it) ); }
429 pybind11::capsule pyCapsuleGrpEl ( npVals, [](
void *f ) {
float* foo =
reinterpret_cast< float*
> ( f );
delete foo; } );
432 pybind11::array_t < float > retArr = pybind11::array_t<float> ( { 3, 3 },
433 { 3 *
sizeof(float),
sizeof(
float) },
438 retList.append ( retArr );
443 },
"This function returns the group elements as rotation matrices of any point group described by the detected axes.", pybind11::arg (
"axList" ), pybind11::arg (
"groupType" ) =
"", pybind11::arg(
"matrixTolerance" ) = 0.05 )
445 .def (
"getMapCOMProcessChange",
449 std::vector< proshade_double > vals =
self.getMapCOMProcessChange ();
452 float* npVals =
new float[
static_cast<unsigned int> ( 3 )];
456 for ( proshade_unsign iter = 0; iter < 3; iter++ ) { npVals[iter] =
static_cast< float > ( vals.at(iter) ); }
459 pybind11::capsule pyCapsuleSymShiftDat ( npVals, [](
void *f ) {
float* foo =
reinterpret_cast< float*
> ( f );
delete foo; } );
462 pybind11::array_t < float > retArr = pybind11::array_t<float> ( {
static_cast<int> ( vals.size() ) },
465 pyCapsuleSymShiftDat );
469 },
"This function returns the shift in Angstrom applied to the internal map representation in order to align its COM with the centre of box." )
473 .def (
"getBestRotationMapPeaksEulerAngles",
477 std::vector< proshade_double > vals =
self.getBestRotationMapPeaksEulerAngles ( settings );
480 float* npVals =
new float[
static_cast<unsigned int> (vals.size())];
484 for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( vals.size() ); iter++ ) { npVals[iter] =
static_cast< float > ( vals.at(iter) ); }
487 pybind11::capsule pyCapsuleEuPeak ( npVals, [](
void *f ) {
float* foo =
reinterpret_cast< float*
> ( f );
delete foo; } );
490 pybind11::array_t < float > retArr = pybind11::array_t<float> ( {
static_cast<unsigned int> (vals.size()) },
497 },
"This function returns a vector of three floats, the three Euler angles of the best peak in the rotation map.", pybind11::arg (
"settings" ) )
498 .def (
"getBestRotationMapPeaksRotationMatrix",
502 std::vector< proshade_double > vals =
self.getBestRotationMapPeaksEulerAngles ( settings );
505 proshade_double* retMat =
new proshade_double[9];
510 pybind11::capsule pyCapsuleRMPeak ( retMat, [](
void *f ) { proshade_double* foo =
reinterpret_cast< proshade_double*
> ( f );
delete foo; } );
513 pybind11::array_t < proshade_double > retArr = pybind11::array_t<proshade_double> ( { 3, 3 },
514 { 3 *
sizeof(proshade_double),
sizeof(proshade_double) },
520 },
"This function returns a rotation matrix representing the best peak in the rotation map.", pybind11::arg (
"settings" ) )
525 .def (
"getOverlayTranslations",
529 std::vector< proshade_double > vals =
self.getBestTranslationMapPeaksAngstrom ( staticStructure );
532 pybind11::dict retDict;
533 pybind11::list rotCen, toOverlay;
536 rotCen.append (
self.originalPdbRotCenX );
537 rotCen.append (
self.originalPdbRotCenY );
538 rotCen.append (
self.originalPdbRotCenZ );
540 toOverlay.append (
self.originalPdbTransX );
541 toOverlay.append (
self.originalPdbTransY );
542 toOverlay.append (
self.originalPdbTransZ );
545 retDict[ pybind11::handle ( pybind11::str (
"centreOfRotation" ).ptr ( ) ) ] = rotCen;
546 retDict[ pybind11::handle ( pybind11::str (
"rotCenToOverlay" ).ptr ( ) ) ] = toOverlay;
550 },
"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" ) )
551 .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 (
"trsX" ), pybind11::arg (
"trsY" ), pybind11::arg (
"trsZ" ) )
554 .def (
"findSHIndex",
558 proshade_signed index = seanindex (
static_cast< int > ( order ),
559 static_cast< int > ( band ),
560 static_cast< int > (
self.spheres[shell]->getLocalBandwidth() ) );
564 },
"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." )
565 .def (
"getSphericalHarmonics",
569 std::complex<proshade_double>* npVals =
new std::complex<proshade_double>[
static_cast<proshade_unsign
> (
self.noSpheres * pow(
self.maxShellBand, 2.0 ) )];
573 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 ); }
576 proshade_signed pyPosSH;
577 proshade_signed pyPos;
580 for ( proshade_signed shIt = 0; shIt < static_cast<proshade_signed> (
self.noSpheres ); shIt++ )
582 for ( proshade_signed bnd = 0; bnd < static_cast<proshade_signed> (
self.spheres[shIt]->getLocalBandwidth() ); bnd++ )
584 for ( proshade_signed order = -bnd; order <= bnd; order++ )
586 pyPosSH = (
static_cast< proshade_signed
> ( shIt ) *
static_cast< proshade_signed
> ( std::pow (
self.maxShellBand, 2 ) ) );
587 pyPos = seanindex (
static_cast< int > ( order ),
588 static_cast< int > ( bnd ),
589 static_cast< int > (
self.spheres[shIt]->getLocalBandwidth() ) );
590 npVals[pyPosSH+pyPos].real (
self.sphericalHarmonics[shIt][pyPos][0] );
591 npVals[pyPosSH+pyPos].imag (
self.sphericalHarmonics[shIt][pyPos][1] );
597 pybind11::capsule pyCapsuleSHs ( npVals, [](
void *f ) { std::complex<proshade_double>* foo =
reinterpret_cast< std::complex<proshade_double>*
> ( f );
delete foo; } );
600 pybind11::array_t < std::complex<proshade_double> > retArr = pybind11::array_t < std::complex<proshade_double > > (
601 {
static_cast<int> (
self.noSpheres ),
static_cast<int> ( pow (
self.maxShellBand, 2.0 ) ) },
602 {
sizeof ( std::complex < proshade_double > ) *
static_cast< proshade_unsign
> ( std::pow (
static_cast< proshade_unsign
> (
self.maxShellBand ), 2 ) ),
sizeof ( std::complex < proshade_double > ) },
608 },
"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 maxEMatDim**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." )
613 std::complex<proshade_double>* npVals =
new std::complex < proshade_double >[
static_cast<proshade_unsign
> (
self.maxEMatDim * ( (
self.maxEMatDim * 2 ) + 1 ) * ( (
self.maxEMatDim * 2 ) + 1 ) )];
617 proshade_signed index = 0;
620 for ( proshade_unsign iter = 0; iter < static_cast < proshade_unsign > (
self.maxEMatDim * ( (
self.maxEMatDim * 2 ) + 1 ) * ( (
self.maxEMatDim * 2 ) + 1 ) ); iter++ ) { npVals[iter].real ( 0.0 ); npVals[iter].imag ( 0.0 ); }
623 proshade_double emReal, emImag;
624 for ( proshade_signed bandIter = 0; bandIter < static_cast< proshade_signed > (
self.maxEMatDim ); bandIter++ )
626 for ( proshade_signed order1 = 0; order1 < ( ( bandIter * 2 ) + 1 ); order1++ )
628 for ( proshade_signed order2 = 0; order2 < ( ( bandIter * 2 ) + 1 ); order2++ )
630 index = order2 + ( (
static_cast< proshade_signed
> (
self.maxEMatDim ) * 2 ) + 1 ) * ( order1 + ( (
static_cast< proshade_signed
> (
self.maxEMatDim ) * 2 ) + 1 ) * bandIter );
631 self.getEMatrixValue(
static_cast< proshade_unsign
> ( bandIter ),
static_cast< proshade_unsign
> ( order1 ),
static_cast< proshade_unsign
> ( order2 ), &emReal, &emImag );
632 npVals[index].real ( emReal );
633 npVals[index].imag ( emImag );
639 pybind11::capsule pyCapsuleEMs ( npVals, [](
void *f ) { std::complex<proshade_double>* foo =
reinterpret_cast< std::complex<proshade_double>*
> ( f );
delete foo; } );
642 proshade_unsign orderSize = ( (
static_cast< proshade_unsign
> (
self.maxEMatDim ) * 2 ) + 1 );
643 pybind11::array_t < std::complex<proshade_double> > retArr = pybind11::array_t < std::complex<proshade_double > > (
644 {
self.maxEMatDim, orderSize, orderSize },
645 {
sizeof ( std::complex < proshade_double > ) * orderSize * orderSize,
646 sizeof ( std::complex < proshade_double > ) * orderSize,
647 sizeof ( std::complex < proshade_double > ) },
653 },
"This function returns the weighted 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 work, the order indexing starts from 0 - i.e. array[1][0][0] means band 1 ; order1 = -1 and order2 = -1." )
654 .def (
"getSO3Coefficients",
658 std::complex<proshade_double>* npVals =
new std::complex < proshade_double >[
static_cast<proshade_unsign
> (
self.maxEMatDim * ( (
self.maxEMatDim * 2 ) + 1 ) * ( (
self.maxEMatDim * 2 ) + 1 ) )];
662 proshade_signed index = 0;
665 for ( proshade_unsign iter = 0; iter < static_cast < proshade_unsign > (
self.maxEMatDim * ( (
self.maxEMatDim * 2 ) + 1 ) * ( (
self.maxEMatDim * 2 ) + 1 ) ); iter++ ) { npVals[iter].real ( 0.0 ); npVals[iter].imag ( 0.0 ); }
668 for ( proshade_signed bandIter = 0; bandIter < static_cast< proshade_signed > (
self.maxEMatDim ); bandIter++ )
670 for ( proshade_signed order1 = 0; order1 < ( ( bandIter * 2 ) + 1 ); order1++ )
672 for ( proshade_signed order2 = 0; order2 < ( ( bandIter * 2 ) + 1 ); order2++ )
674 index = order2 + ( (
static_cast< proshade_signed
> (
self.maxEMatDim ) * 2 ) + 1 ) * ( order1 + ( (
static_cast< proshade_signed
> (
self.maxEMatDim ) * 2 ) + 1 ) * bandIter );
675 npVals[index].real (
self.so3Coeffs[
self.so3CoeffsArrayIndex (
static_cast< proshade_signed
> ( order1 - bandIter ), order2 - bandIter, bandIter )][0] );
676 npVals[index].imag (
self.so3Coeffs[
self.so3CoeffsArrayIndex (
static_cast< proshade_signed
> ( order1 - bandIter ), order2 - bandIter, bandIter )][1] );
682 pybind11::capsule pyCapsuleSOCoeffs ( npVals, [](
void *f ) { std::complex<proshade_double>* foo =
reinterpret_cast< std::complex<proshade_double>*
> ( f );
delete foo; } );
685 proshade_unsign orderSize = ( (
static_cast< proshade_unsign
> (
self.maxEMatDim ) * 2 ) + 1 );
686 pybind11::array_t < std::complex<proshade_double> > retArr = pybind11::array_t < std::complex<proshade_double > > (
687 {
self.maxEMatDim, orderSize, orderSize },
688 {
sizeof ( std::complex < proshade_double > ) * orderSize * orderSize,
689 sizeof ( std::complex < proshade_double > ) * orderSize,
690 sizeof ( std::complex < proshade_double > ) },
696 },
"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." )
697 .def (
"getRotationFunctionMap",
701 std::complex<proshade_double>* npVals =
new std::complex < proshade_double >[
static_cast<proshade_unsign
> ( (
self.maxEMatDim * 2 ) * (
self.maxEMatDim * 2 ) * (
self.maxEMatDim * 2 ) )];
705 for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( (
self.maxEMatDim * 2 ) * (
self.maxEMatDim * 2 ) * (
self.maxEMatDim * 2 ) ); iter++ ) { npVals[iter].real (
self.so3CoeffsInverse[iter][0] ); npVals[iter].imag (
self.so3CoeffsInverse[iter][1] ); }
708 pybind11::capsule pyCapsuleRotMap ( npVals, [](
void *f ) { std::complex<proshade_double>* foo =
reinterpret_cast< std::complex<proshade_double>*
> ( f );
delete foo; } );
711 pybind11::array_t < std::complex < proshade_double > > retArr = pybind11::array_t < std::complex < proshade_double > > (
712 { (
self.maxEMatDim * 2 ), (
self.maxEMatDim * 2 ), (
self.maxEMatDim * 2 ) },
713 { (
self.maxEMatDim * 2 ) * (
self.maxEMatDim * 2 ) *
sizeof(std::complex < proshade_double >),
714 (
self.maxEMatDim * 2 ) *
sizeof(std::complex < proshade_double >),
715 sizeof(std::complex < proshade_double >) },
721 },
"This function returns the (self) rotation function as a three-dimensional map of complex numbers." )
722 .def (
"getRotationMatrixFromSOFTCoordinates",
726 proshade_double* npVals =
new proshade_double[9];
730 proshade_double eulA, eulB, eulG;
739 pybind11::capsule pyCapsuleRMSoft ( npVals, [](
void *f ) { proshade_double* foo =
reinterpret_cast< proshade_double*
> ( f );
delete foo; } );
742 pybind11::array_t < proshade_double > retArr = pybind11::array_t < proshade_double > ( { 3, 3 },
743 { 3 *
sizeof(proshade_double),
sizeof(proshade_double) },
749 },
"This function converts a given rotation function map position onto the corresponding rotation matrix.", pybind11::arg (
"xPos" ), pybind11::arg (
"yPos" ), pybind11::arg (
"zPos" ) )
750 .def (
"getTranslationFunctionMap",
754 std::complex<proshade_double>* npVals =
new std::complex < proshade_double >[
static_cast<proshade_unsign
> (
self.getXDim() *
self.getYDim() *
self.getZDim() )];
758 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] ); }
761 pybind11::capsule pyCapsuleTrsMap ( npVals, [](
void *f ) { std::complex<proshade_double>* foo =
reinterpret_cast< std::complex<proshade_double>*
> ( f );
delete foo; } );
764 pybind11::array_t < std::complex < proshade_double > > retArr = pybind11::array_t < std::complex < proshade_double > > (
765 {
self.getXDim(),
self.getYDim(),
self.getZDim() },
766 {
self.getYDim() *
self.getZDim() *
sizeof(std::complex < proshade_double >),
767 self.getZDim() *
sizeof(std::complex < proshade_double >),
768 sizeof(std::complex < proshade_double >) },
774 },
"This function returns the translation function as a three-dimensional map of complex numbers." )
813 .def (
"__repr__", [] ( ) {
return "<ProSHADE_data class object> (This class contains all information, results and available functionalities for a structure)"; } );
816 pyProSHADE.def (
"computeGroupElementsForGroup",
817 [] ( proshade_double x, proshade_double y, proshade_double z, proshade_unsign fold ) -> pybind11::array_t < proshade_double >
823 proshade_double* npVals =
new proshade_double[
static_cast<proshade_unsign
> ( retVec.size() ) * 9];
827 for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( retVec.size() ); iter++ ) {
for ( proshade_unsign it = 0; it < 9; it++ ) { npVals[(iter*9)+it] = retVec.at(iter).at(it); } }
830 pybind11::capsule pyCapsuleGrEls ( npVals, [](
void *f ) { proshade_double* foo =
reinterpret_cast< proshade_double*
> ( f );
delete foo; } );
833 pybind11::array_t < proshade_double > retArr = pybind11::array_t < proshade_double > (
834 {
static_cast< int > ( retVec.size() ), 3, 3 },
835 { 9 *
sizeof(proshade_double), 3 *
sizeof(proshade_double),
sizeof(proshade_double) },
841 },
"This function takes the axis and fold and computes all resulting group elements (as rotation matrices), returning them in numpy.ndarray.", pybind11::arg (
"x" ), pybind11::arg (
"y" ), pybind11::arg (
"z" ), pybind11::arg (
"fold" ) );
842 pyProSHADE.def (
"joinElementsFromDifferentGroups",
843 [] ( pybind11::array_t < proshade_double, pybind11::array::c_style | pybind11::array::forcecast > first,
844 pybind11::array_t < proshade_double, pybind11::array::c_style | pybind11::array::forcecast > second,
845 proshade_double matrixTolerance,
846 bool combine ) -> pybind11::array_t < proshade_double >
849 pybind11::buffer_info buf1 = first.request();
850 pybind11::buffer_info buf2 = second.request();
853 if ( buf1.ndim != 3 || buf2.ndim != 3 )
855 std::cerr <<
"Function joinElementsFromDifferentGroups() arguments first and second should be numpy.ndarrays with 3 dimensions indexed as follos: first[elementNumber][elementRotationMatrixRow][elementRotationMatrixColumn] - the same format as returned by the computeGroupElementsForGroup() function." << std::endl;
856 return ( pybind11::array_t < proshade_double > () );
860 std::vector< std::vector< proshade_double > > fVec, sVec;
862 proshade_double* dataPtr1 =
reinterpret_cast < proshade_double*
> ( buf1.ptr );
863 for ( proshade_unsign elIt = 0; elIt < static_cast<proshade_unsign> ( buf1.shape.at(0) ); elIt++ )
865 std::vector< proshade_double > rotMat;
866 for ( proshade_unsign rowIt = 0; rowIt < static_cast<proshade_unsign> ( buf1.shape.at(1) ); rowIt++ )
868 for ( proshade_unsign colIt = 0; colIt < static_cast<proshade_unsign> ( buf1.shape.at(2) ); colIt++ )
876 proshade_double* dataPtr2 =
reinterpret_cast < proshade_double*
> ( buf2.ptr );
877 for ( proshade_unsign elIt = 0; elIt < static_cast<proshade_unsign> ( buf2.shape.at(0) ); elIt++ )
879 std::vector< proshade_double > rotMat;
880 for ( proshade_unsign rowIt = 0; rowIt < static_cast<proshade_unsign> ( buf2.shape.at(1) ); rowIt++ )
882 for ( proshade_unsign colIt = 0; colIt < static_cast<proshade_unsign> ( buf2.shape.at(2) ); colIt++ )
895 proshade_double* npVals =
new proshade_double[
static_cast<proshade_unsign
> ( retVec.size() ) * 9];
899 for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( retVec.size() ); iter++ ) {
for ( proshade_unsign it = 0; it < 9; it++ ) { npVals[(iter*9)+it] = retVec.at(iter).at(it); } }
902 pybind11::capsule pyCapsuleGrElsCombo ( npVals, [](
void *f ) { proshade_double* foo =
reinterpret_cast< proshade_double*
> ( f );
delete foo; } );
905 pybind11::array_t < proshade_double > retArr = pybind11::array_t < proshade_double > (
906 {
static_cast< int > ( retVec.size() ), 3, 3 },
907 { 9 *
sizeof(proshade_double), 3 *
sizeof(proshade_double),
sizeof(proshade_double) },
909 pyCapsuleGrElsCombo );
913 },
"This function takes the axis and fold and computes all resulting group elements (as rotation matrices), returning them in numpy.ndarray.", pybind11::arg (
"first" ), pybind11::arg (
"second" ), pybind11::arg (
"matrixTolerance" ), pybind11::arg (
"combine" ) );
914 pyProSHADE.def (
"getRotationFunctionSpheres",
921 std::cerr <<
"Function getRotationFunctionSpheres() was called before rotation function was computed. Please compute the rotation function by calling the computeRotationFunction() function first. Returning empty array." << std::endl;
922 return ( pybind11::array_t < proshade_double > () );
928 std::cerr <<
"The fold for which the rotation function should be mapped is not in the allowed range (i.e. 2+). Returning empty array." << std::endl;
929 return ( pybind11::array_t < proshade_double > () );
933 std::vector<ProSHADE_internal_spheres::ProSHADE_rotFun_sphere*> sphereMappedRotFun;
936 proshade_double soughtAngle = 0.0;
937 for ( proshade_double angIt = 1.0; angIt < static_cast < proshade_double > ( fold ); angIt += 1.0 )
940 soughtAngle = angIt * ( 2.0 * M_PI /
static_cast<proshade_double
> ( fold ) );
944 M_PI /
static_cast< proshade_double
> ( dataObj->
getMaxBand() ),
948 static_cast<proshade_unsign
> ( angIt - 1.0 ) ) );
951 sphereMappedRotFun.at(
static_cast < size_t > ( angIt - 1.0 ))->interpolateSphereValues ( dataObj->
getInvSO3Coeffs ( ) );
955 proshade_double* npVals =
new proshade_double[(dataObj->
getMaxBand() * 2) * (dataObj->
getMaxBand() * 2) * (
static_cast< proshade_unsign
> ( fold ) - 1)];
959 std::vector< std::vector < proshade_double > > vls;
961 for (
size_t sphIt = 0; sphIt < sphereMappedRotFun.size(); sphIt++ )
964 vls = sphereMappedRotFun.at(sphIt)->getCopyOfValues ();
966 for (
size_t lonIt = 0; lonIt < (dataObj->
getMaxBand() * 2); lonIt++ )
968 for (
size_t latIt = 0; latIt < (dataObj->
getMaxBand() * 2); latIt++ )
970 valIter = latIt +
static_cast< size_t > ( dataObj->
getMaxBand() * 2 ) * ( lonIt +
static_cast< size_t > ( dataObj->
getMaxBand() * 2 ) * sphIt );
971 npVals[valIter] = vls.at(lonIt).at(latIt);
977 pybind11::capsule pyCapsuleRFMapSph ( npVals, [](
void *f ) { proshade_double* foo =
reinterpret_cast< proshade_double*
> ( f );
delete foo; } );
980 pybind11::array_t < proshade_double > retArr = pybind11::array_t < proshade_double > (
981 {
static_cast< int > ( sphereMappedRotFun.size() ),
static_cast< int > (dataObj->
getMaxBand() * 2),
static_cast< int > (dataObj->
getMaxBand() * 2) },
982 { (dataObj->
getMaxBand() * 2) * (dataObj->
getMaxBand() * 2) *
sizeof(proshade_double), (dataObj->
getMaxBand() * 2) *
sizeof(proshade_double),
sizeof(proshade_double) },
988 },
"This function allows access to the (self-)rotation function mapping to spheres with the angle being the radius and the longitude/latitude angles being the axis of the rotation for a particular fold (i.e. for all angles required by the fold except for the angle 0)." );