nmie-pybind11.cc 8.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180
  1. //**********************************************************************************//
  2. // Copyright (C) 2009-2018 Ovidio Pena <ovidio@bytesfall.com> //
  3. // Copyright (C) 2013-2018 Konstantin Ladutenko <kostyfisik@gmail.com> //
  4. // //
  5. // This file is part of scattnlay //
  6. // //
  7. // This program is free software: you can redistribute it and/or modify //
  8. // it under the terms of the GNU General Public License as published by //
  9. // the Free Software Foundation, either version 3 of the License, or //
  10. // (at your option) any later version. //
  11. // //
  12. // This program is distributed in the hope that it will be useful, //
  13. // but WITHOUT ANY WARRANTY; without even the implied warranty of //
  14. // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the //
  15. // GNU General Public License for more details. //
  16. // //
  17. // The only additional remark is that we expect that all publications //
  18. // describing work using this software, or all commercial products //
  19. // using it, cite at least one of the following references: //
  20. // [1] O. Pena and U. Pal, "Scattering of electromagnetic radiation by //
  21. // a multilayered sphere," Computer Physics Communications, //
  22. // vol. 180, Nov. 2009, pp. 2348-2354. //
  23. // [2] K. Ladutenko, U. Pal, A. Rivera, and O. Pena-Rodriguez, "Mie //
  24. // calculation of electromagnetic near-field for a multilayered //
  25. // sphere," Computer Physics Communications, vol. 214, May 2017, //
  26. // pp. 225-230. //
  27. // //
  28. // You should have received a copy of the GNU General Public License //
  29. // along with this program. If not, see <http://www.gnu.org/licenses/>. //
  30. //**********************************************************************************//
  31. #include "nmie.hpp"
  32. #include "nmie-impl.hpp"
  33. #include "nmie-precision.hpp"
  34. #include <array>
  35. #include <tuple>
  36. #include <algorithm>
  37. #include <cstdio>
  38. #include <cstdlib>
  39. #include <stdexcept>
  40. #include <iostream>
  41. #include <iomanip>
  42. #include <vector>
  43. #include <pybind11/pybind11.h>
  44. #include <pybind11/numpy.h>
  45. #include <pybind11/complex.h>
  46. namespace py = pybind11;
  47. py::array_t< std::complex<double>> VectorComplex2Py(const std::vector<std::complex<double> > &c_x) {
  48. auto py_x = py::array_t< std::complex<double>>(c_x.size());
  49. auto py_x_buffer = py_x.request();
  50. std::complex<double> *py_x_ptr = ( std::complex<double> *) py_x_buffer.ptr;
  51. std::memcpy(py_x_ptr,c_x.data(),c_x.size()*sizeof( std::complex<double>));
  52. return py_x;
  53. }
  54. // https://stackoverflow.com/questions/17294629/merging-flattening-sub-vectors-into-a-single-vector-c-converting-2d-to-1d
  55. template <typename T>
  56. std::vector<T> flatten(const std::vector<std::vector<T>>& v) {
  57. std::size_t total_size = 0;
  58. for (const auto& sub : v)
  59. total_size += sub.size(); // I wish there was a transform_accumulate
  60. std::vector<T> result;
  61. result.reserve(total_size);
  62. for (const auto& sub : v)
  63. result.insert(result.end(), sub.begin(), sub.end());
  64. return result;
  65. }
  66. template <typename T>
  67. py::array VectorVector2Py(const std::vector<std::vector<T > > &x) {
  68. size_t dim1 = x.size();
  69. size_t dim2 = x[0].size();
  70. auto result = flatten(x);
  71. // https://github.com/tdegeus/pybind11_examples/blob/master/04_numpy-2D_cpp-vector/example.cpp
  72. size_t ndim = 2;
  73. std::vector<size_t> shape = { dim1 , dim2 };
  74. std::vector<size_t> strides = { sizeof(T)*dim2 , sizeof(T) };
  75. // return 2-D NumPy array
  76. return py::array(py::buffer_info(
  77. result.data(), /* data as contiguous array */
  78. sizeof(T), /* size of one scalar */
  79. py::format_descriptor<T>::format(), /* data type */
  80. ndim, /* number of dimensions */
  81. shape, /* shape of the matrix */
  82. strides /* strides for each axis */
  83. ));
  84. }
  85. template <typename T>
  86. std::vector<T> Py2Vector(const py::array_t<T> &py_x) {
  87. std::vector<T> c_x(py_x.size());
  88. std::memcpy(c_x.data(), py_x.data(), py_x.size()*sizeof(T));
  89. return c_x;
  90. }
  91. py::tuple py_ScattCoeffs(const py::array_t<double, py::array::c_style | py::array::forcecast> &py_x,
  92. const py::array_t< std::complex<double>, py::array::c_style | py::array::forcecast> &py_m,
  93. const int nmax=-1, const int pl=-1) {
  94. auto c_x = Py2Vector<double>(py_x);
  95. auto c_m = Py2Vector< std::complex<double> >(py_m);
  96. int terms = 0;
  97. std::vector<std::complex<double> > c_an, c_bn;
  98. int L = py_x.size();
  99. terms = nmie::ScattCoeffs( L, pl, c_x, c_m, nmax, c_an, c_bn);
  100. return py::make_tuple(terms, VectorComplex2Py(c_an), VectorComplex2Py(c_bn));
  101. }
  102. py::tuple py_scattnlay(const py::array_t<double, py::array::c_style | py::array::forcecast> &py_x,
  103. const py::array_t< std::complex<double>, py::array::c_style | py::array::forcecast> &py_m,
  104. const py::array_t<double, py::array::c_style | py::array::forcecast> &py_theta,
  105. const int nmax=-1, const int pl=-1) {
  106. auto c_x = Py2Vector<double>(py_x);
  107. auto c_m = Py2Vector< std::complex<double> >(py_m);
  108. auto c_theta = Py2Vector<double>(py_theta);
  109. int L = py_x.size(), nTheta = c_theta.size(), terms;
  110. double Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo;
  111. std::vector<std::complex<double> > c_S1, c_S2;
  112. terms = nmie::nMie(L, pl, c_x, c_m, nTheta, c_theta, nmax, &Qext, &Qsca, &Qabs, &Qbk, &Qpr, &g, &Albedo, c_S1, c_S2);
  113. return py::make_tuple(terms, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo,
  114. VectorComplex2Py(c_S1), VectorComplex2Py(c_S2));
  115. }
  116. py::tuple py_fieldnlay(const py::array_t<double, py::array::c_style | py::array::forcecast> &py_x,
  117. const py::array_t< std::complex<double>, py::array::c_style | py::array::forcecast> &py_m,
  118. const py::array_t<double, py::array::c_style | py::array::forcecast> &py_Xp,
  119. const py::array_t<double, py::array::c_style | py::array::forcecast> &py_Yp,
  120. const py::array_t<double, py::array::c_style | py::array::forcecast> &py_Zp,
  121. const int nmax=-1, const int pl=-1) {
  122. auto c_x = Py2Vector<double>(py_x);
  123. auto c_m = Py2Vector< std::complex<double> >(py_m);
  124. auto c_Xp = Py2Vector<double>(py_Xp);
  125. auto c_Yp = Py2Vector<double>(py_Yp);
  126. auto c_Zp = Py2Vector<double>(py_Zp);
  127. unsigned int ncoord = py_Xp.size();
  128. std::vector<std::vector<std::complex<double> > > E(ncoord);
  129. std::vector<std::vector<std::complex<double> > > H(ncoord);
  130. for (auto& f : E) f.resize(3);
  131. for (auto& f : H) f.resize(3);
  132. int L = py_x.size(), terms;
  133. terms = nmie::nField(L, pl, c_x, c_m, nmax, ncoord, c_Xp, c_Yp, c_Zp, E, H);
  134. auto py_E = VectorVector2Py<std::complex<double> >(E);
  135. auto py_H = VectorVector2Py<std::complex<double> >(H);
  136. return py::make_tuple(terms, py_E, py_H);
  137. }
  138. PYBIND11_MODULE(example, m) {
  139. m.doc() = "pybind11 example plugin"; // optional module docstring
  140. m.def("scattcoeffs_", &py_ScattCoeffs, "test",
  141. py::arg("x"), py::arg("m"),
  142. py::arg("nmax")=-1, py::arg("pl")=-1);
  143. m.def("scattnlay_", &py_scattnlay, "test",
  144. py::arg("x"), py::arg("m"),
  145. py::arg("theta")=py::array_t<double>(0),
  146. py::arg("nmax")=-1, py::arg("pl")=-1);
  147. m.def("fieldnlay_", &py_fieldnlay, "test",
  148. py::arg("x"), py::arg("m"),
  149. py::arg("xp"),py::arg("yp"),py::arg("zp"),
  150. py::arg("nmax")=-1, py::arg("pl")=-1);
  151. }