#include #include "nmie.hpp" #include "nmie-basic.hpp" #include "nmie-nearfield.hpp" #include #include #include namespace py = pybind11; template std::vector Py2Vector(const py::array_t &py_x) { std::vector c_x(py_x.size()); std::memcpy(c_x.data(), py_x.data(), py_x.size()*sizeof(T)); return c_x; } // https://github.com/pybind/pybind11/issues/1042#issuecomment-508582847 //template //inline py::array_t Vector2Py(Sequence&& seq) { // // Move entire object to heap (Ensure is moveable!). Memory handled via Python capsule // Sequence* seq_ptr = new Sequence(std::move(seq)); // auto capsule = py::capsule(seq_ptr, [](void* p) { delete reinterpret_cast(p); }); // return py::array(seq_ptr->size(), // shape of array // seq_ptr->data(), // c-style contiguous strides for Sequence // capsule // numpy array references this parent // ); //} template inline py::array_t Vector2Py(const std::vector& seq) { return py::array(seq.size(), seq.data()); } template py::array_t< std::complex> VectorComplex2Py(const std::vector > &cf_x) { auto c_x = nmie::ConvertComplexVector(cf_x); auto py_x = py::array_t< std::complex>(c_x.size()); auto py_x_buffer = py_x.request(); auto *py_x_ptr = (std::complex *) py_x_buffer.ptr; std::memcpy(py_x_ptr, c_x.data(), c_x.size()*sizeof(std::complex)); return py_x; } // https://stackoverflow.com/questions/17294629/merging-flattening-sub-vectors-into-a-single-vector-c-converting-2d-to-1d template std::vector flatten(const std::vector> &v) { std::size_t total_size = 0; for (const auto &sub : v) total_size += sub.size(); // I wish there was a transform_accumulate std::vector result; result.reserve(total_size); for (const auto &sub : v) result.insert(result.end(), sub.begin(), sub.end()); return result; } template py::array Vector2DComplex2Py(const std::vector > &x) { size_t dim1 = x.size(); size_t dim2 = x[0].size(); auto result = flatten(x); // https://github.com/tdegeus/pybind11_examples/blob/master/04_numpy-2D_cpp-vector/example.cpp size_t ndim = 2; std::vector shape = {dim1, dim2}; std::vector strides = {sizeof(T)*dim2, sizeof(T)}; // return 2-D NumPy array return py::array(py::buffer_info( result.data(), /* data as contiguous array */ sizeof(T), /* size of one scalar */ py::format_descriptor::format(), /* data type */ ndim, /* number of dimensions */ shape, /* shape of the matrix */ strides /* strides for each axis */ )); } namespace nmie { template class PyMultiLayerMie : public MultiLayerMie { public: template py::array_t> GetS1(); template py::array_t> GetS2(); template py::array_t> GetAn(); template py::array_t> GetBn(); template py::array GetLayerAn(); template py::array GetLayerBn(); template py::array GetLayerCn(); template py::array GetLayerDn(); template py::array GetFieldE(); template py::array GetFieldH(); template py::array_t GetFieldEabs(); template py::array_t GetFieldHabs(); template void SetLayersSize(const py::array_t &py_layer_size); template void SetLayersIndex(const py::array_t, py::array::c_style | py::array::forcecast> &py_index); template void SetAngles(const py::array_t &py_angles); void SetFieldCoords(const py::array_t &py_Xp, const py::array_t &py_Yp, const py::array_t &py_Zp); }; // Python interface template template void PyMultiLayerMie::SetLayersSize( const py::array_t &py_layer_size) { auto layer_size_dp = Py2Vector(py_layer_size); this->MultiLayerMie::SetLayersSize(ConvertVector(layer_size_dp)); } template template void PyMultiLayerMie::SetLayersIndex( const py::array_t, py::array::c_style | py::array::forcecast> &py_index) { auto index_dp = Py2Vector >(py_index); this->MultiLayerMie::SetLayersIndex(ConvertComplexVector(index_dp)); } template template void PyMultiLayerMie::SetAngles(const py::array_t &py_angles) { auto angles_dp = Py2Vector(py_angles); this->MultiLayerMie::SetAngles(ConvertVector(angles_dp)); } template template py::array_t> PyMultiLayerMie::GetS1() { return VectorComplex2Py( this->MultiLayerMie::GetS1() ); } template template py::array_t> PyMultiLayerMie::GetS2() { return VectorComplex2Py( this->MultiLayerMie::GetS2() ); } template template py::array_t PyMultiLayerMie::GetFieldEabs() { return Vector2Py(ConvertVector( this->MultiLayerMie::GetFieldEabs() )); } template template py::array_t PyMultiLayerMie::GetFieldHabs() { return Vector2Py(ConvertVector( this->MultiLayerMie::GetFieldHabs() )); } template template py::array_t> PyMultiLayerMie::GetAn() { return VectorComplex2Py( this->MultiLayerMie::GetAn() ); } template template py::array_t> PyMultiLayerMie::GetBn() { return VectorComplex2Py( this->MultiLayerMie::GetBn() ); } template template py::array PyMultiLayerMie::GetFieldE() { return Vector2DComplex2Py >( ConvertComplexVectorVector( this->MultiLayerMie::GetFieldE() ) ); } template template py::array PyMultiLayerMie::GetFieldH() { return Vector2DComplex2Py >( ConvertComplexVectorVector( this->MultiLayerMie::GetFieldH() ) ); } template template py::array PyMultiLayerMie::GetLayerAn() { return Vector2DComplex2Py >( ConvertComplexVectorVector( this->MultiLayerMie::GetLayerAn() ) ); } template template py::array PyMultiLayerMie::GetLayerBn() { return Vector2DComplex2Py >( ConvertComplexVectorVector( this->MultiLayerMie::GetLayerBn() ) ); } template template py::array PyMultiLayerMie::GetLayerCn() { return Vector2DComplex2Py >( ConvertComplexVectorVector( this->MultiLayerMie::GetLayerCn() ) ); } template template py::array PyMultiLayerMie::GetLayerDn() { return Vector2DComplex2Py >( ConvertComplexVectorVector( this->MultiLayerMie::GetLayerDn() ) ); } template void PyMultiLayerMie::SetFieldCoords( const py::array_t &py_Xp, const py::array_t &py_Yp, const py::array_t &py_Zp) { auto c_Xp = Py2Vector(py_Xp); auto c_Yp = Py2Vector(py_Yp); auto c_Zp = Py2Vector(py_Zp); this->MultiLayerMie::SetFieldCoords({ConvertVector(c_Xp), ConvertVector(c_Yp), ConvertVector(c_Zp)}); } } template void declare_nmie(py::module &m, const std::string &typestr) { using mie_typed = nmie::PyMultiLayerMie; std::string pyclass_name = std::string("mie") + typestr; py::class_(m, pyclass_name.c_str(), py::buffer_protocol(), py::dynamic_attr()) .def(py::init<>()) .def("GetPECLayer", &mie_typed::GetPECLayer) .def("SetPECLayer", &mie_typed::SetPECLayer) .def("SetMaxTerms", &mie_typed::SetMaxTerms) .def("GetMaxTerms", &mie_typed::GetMaxTerms) .def("SetModeNmaxAndType", &mie_typed::SetModeNmaxAndType) .def("RunMieCalculation", &mie_typed::RunMieCalculation) .def("calcScattCoeffs", &mie_typed::calcScattCoeffs) .def("calcExpanCoeffs", &mie_typed::calcExpanCoeffs) .def("RunFieldCalculation", &mie_typed::RunFieldCalculation) .def("RunFieldCalculationPolar", &mie_typed::RunFieldCalculationPolar) .def("GetPECLayer", &mie_typed::GetPECLayer) .def("SetLayersSize", static_cast< void (mie_typed::*) (const py::array_t&)> (&mie_typed::SetLayersSize) ) .def("SetLayersIndex", static_cast< void (mie_typed::*) (const py::array_t, py::array::c_style | py::array::forcecast> &)> (&mie_typed::SetLayersIndex) ) .def("SetAngles", static_cast< void (mie_typed::*) (const py::array_t&)> (&mie_typed::SetAngles) ) .def("SetFieldCoords", static_cast< void (mie_typed::*) ( const py::array_t &py_Xp, const py::array_t &py_Yp, const py::array_t &py_Zp )> (&mie_typed::SetFieldCoords) ) .def("GetQext", &mie_typed::GetQext) .def("GetQsca", &mie_typed::GetQsca) .def("GetQabs", &mie_typed::GetQabs) .def("GetQpr", &mie_typed::GetQpr) .def("GetQbk", &mie_typed::GetQbk) .def("GetAsymmetryFactor", &mie_typed::GetAsymmetryFactor) .def("GetAlbedo", &mie_typed::GetAlbedo) .def("GetS1", &mie_typed::GetS1) .def("GetS2", &mie_typed::GetS2) .def("GetAn", &mie_typed::GetAn) .def("GetBn", &mie_typed::GetBn) .def("GetFieldE", &mie_typed::GetFieldE) .def("GetFieldH", &mie_typed::GetFieldH) .def("GetFieldEabs", &mie_typed::GetFieldEabs) .def("GetFieldHabs", &mie_typed::GetFieldHabs) .def("GetFieldConvergence", &mie_typed::GetFieldConvergence) .def("GetLayerAn", &mie_typed::GetLayerAn) .def("GetLayerBn", &mie_typed::GetLayerBn) .def("GetLayerCn", &mie_typed::GetLayerCn) .def("GetLayerDn", &mie_typed::GetLayerDn) ; } // wrap as Python module #ifdef MULTI_PRECISION std::string precision_name = "_mp"; PYBIND11_MODULE(scattnlay_mp, m) #else std::string precision_name = "_dp"; PYBIND11_MODULE(scattnlay_dp, m) #endif // MULTI_PRECISION { m.doc() = "The Python version of scattnlay"; declare_nmie(m, precision_name); }