Explorar o código

move to additional file

Konstantin Ladutenko %!s(int64=6) %!d(string=hai) anos
pai
achega
4bce8446d3
Modificáronse 5 ficheiros con 105 adicións e 132 borrados
  1. 5 0
      .gitignore
  2. 1 1
      Makefile
  3. 3 3
      examples/scattcoeffs.py
  4. 96 0
      src/nmie-pybind11.cc
  5. 0 128
      src/nmie.cc

+ 5 - 0
.gitignore

@@ -43,3 +43,8 @@ build*
 tests/python/field.txt
 examples/*.pdf
 examples/*.png
+examples/*.7z
+examples/Egor*
+examples/egor*
+examples/*.txt
+examples/*.pyc

+ 1 - 1
Makefile

@@ -66,7 +66,7 @@ fieldnlay-mp: $(SRCDIR)/nearfield.cc $(SRCDIR)/nmie.cc $(CXX_NMIE_HEADERS)
 
 
 lib:
-	c++ -O3 -Wall -shared -std=c++11 -fPIC `python3 -m pybind11 --includes` $(SRCDIR)/nmie.cc  -lm -o example`python3-config --extension-suffix`
+	c++ -O3 -Wall -shared -std=c++11 -fPIC `python3 -m pybind11 --includes` $(SRCDIR)/nmie.cc $(SRCDIR)/nmie-pybind11.cc  -lm -o example`python3-config --extension-suffix`
 clean:
 	$(PYTHON) setup.py clean
 	$(MAKE) -f $(CURDIR)/debian/rules clean

+ 3 - 3
examples/scattcoeffs.py

@@ -70,8 +70,8 @@ m[:, 3] *= 2.8 + 0.2j
 m[:, 4] *= 1.5 + 0.4j
 
 terms, an, bn = scattcoeffs(x, m, 105)
-terms1, an1, bn1 = example.scattcoeffs(x[0,:], m[0,:])
-print(terms)
+terms1, an1, bn1 = example.scattcoeffs(x[0,:], m[0,:], nmax=10)
+print(an1[:3], bn1[:3])
 print(terms1)
 
 result = np.vstack((x[:, 4], an[:, 0].real, an[:, 0].imag, an[:, 1].real, an[:, 1].imag, an[:, 2].real, an[:, 2].imag,
@@ -96,5 +96,5 @@ try:
     plt.show()
 finally:
     np.savetxt("scattcoeffs.txt", result, fmt = "%.5f")
-    print( result)
+    print( result[0,:])
 

+ 96 - 0
src/nmie-pybind11.cc

@@ -0,0 +1,96 @@
+//**********************************************************************************//
+//    Copyright (C) 2009-2018  Ovidio Pena <ovidio@bytesfall.com>                   //
+//    Copyright (C) 2013-2018  Konstantin Ladutenko <kostyfisik@gmail.com>          //
+//                                                                                  //
+//    This file is part of scattnlay                                                //
+//                                                                                  //
+//    This program is free software: you can redistribute it and/or modify          //
+//    it under the terms of the GNU General Public License as published by          //
+//    the Free Software Foundation, either version 3 of the License, or             //
+//    (at your option) any later version.                                           //
+//                                                                                  //
+//    This program is distributed in the hope that it will be useful,               //
+//    but WITHOUT ANY WARRANTY; without even the implied warranty of                //
+//    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the                 //
+//    GNU General Public License for more details.                                  //
+//                                                                                  //
+//    The only additional remark is that we expect that all publications            //
+//    describing work using this software, or all commercial products               //
+//    using it, cite at least one of the following references:                      //
+//    [1] O. Pena and U. Pal, "Scattering of electromagnetic radiation by           //
+//        a multilayered sphere," Computer Physics Communications,                  //
+//        vol. 180, Nov. 2009, pp. 2348-2354.                                       //
+//    [2] K. Ladutenko, U. Pal, A. Rivera, and O. Pena-Rodriguez, "Mie              //
+//        calculation of electromagnetic near-field for a multilayered              //
+//        sphere," Computer Physics Communications, vol. 214, May 2017,             //
+//        pp. 225-230.                                                              //
+//                                                                                  //
+//    You should have received a copy of the GNU General Public License             //
+//    along with this program.  If not, see <http://www.gnu.org/licenses/>.         //
+//**********************************************************************************//
+#include "nmie.hpp"
+#include "nmie-impl.hpp"
+#include "nmie-precision.hpp"
+#include <array>
+#include <tuple>
+#include <algorithm>
+#include <cstdio>
+#include <cstdlib>
+#include <stdexcept>
+#include <iostream>
+#include <iomanip>
+#include <vector>
+
+#include <pybind11/pybind11.h>
+#include <pybind11/numpy.h>
+#include <pybind11/complex.h>
+
+
+namespace py = pybind11;
+
+
+py::array_t< std::complex<double>> VectorComplex2Py(std::vector<std::complex<double> > c_x) {
+  auto py_x = py::array_t< std::complex<double>>(c_x.size());
+  auto py_x_buffer = py_x.request();
+  std::complex<double> *py_x_ptr    = ( std::complex<double> *) py_x_buffer.ptr;
+  std::memcpy(py_x_ptr,c_x.data(),c_x.size()*sizeof( std::complex<double>));
+  return py_x;
+}
+
+
+std::vector<double> Py2VectorDouble(py::array_t<double> py_x) {
+  std::vector<double> c_x(py_x.size());
+  std::memcpy(c_x.data(), py_x.data(), py_x.size()*sizeof(double));
+  return c_x;
+}
+
+
+std::vector< std::complex<double> > Py2VectorComplex(py::array_t< std::complex<double> > py_x){
+  std::vector< std::complex<double> > c_x(py_x.size());
+  std::memcpy(c_x.data(), py_x.data(), py_x.size()*sizeof( std::complex<double>));
+  return c_x;
+}
+
+
+py::tuple py_ScattCoeffs(py::array_t<double, py::array::c_style | py::array::forcecast> py_x,
+                         py::array_t< std::complex<double>, py::array::c_style | py::array::forcecast> py_m,
+                         const int nmax=-1, const int pl=-1) {
+
+  auto c_x = Py2VectorDouble(py_x);
+  auto c_m = Py2VectorComplex(py_m);
+
+  int terms = 0;
+  std::vector<std::complex<double> > c_an, c_bn;
+  int L = py_x.size();
+  terms = nmie::ScattCoeffs( L, pl, c_x, c_m, nmax, c_an, c_bn);
+  
+  return py::make_tuple(terms, VectorComplex2Py(c_an), VectorComplex2Py(c_bn));
+}
+
+PYBIND11_MODULE(example, m) {
+    m.doc() = "pybind11 example plugin"; // optional module docstring
+
+    m.def("scattcoeffs", &py_ScattCoeffs, "test",
+          py::arg("x"), py::arg("m"),
+          py::arg("nmax")=-1, py::arg("pl")=-1);
+}

+ 0 - 128
src/nmie.cc

@@ -58,17 +58,6 @@
 #include <iomanip>
 #include <vector>
 
-#include <pybind11/pybind11.h>
-#include <pybind11/numpy.h>
-#include <pybind11/complex.h>
-
-
-namespace py = pybind11;
-
-
-
-
-
 namespace nmie {
   
   //**********************************************************************************//
@@ -358,120 +347,3 @@ namespace nmie {
     return 0;
   }
 }  // end of namespace nmie
-
-// // -------------
-// // pure C++ code
-// // -------------
-
-// std::vector<double> length(const std::vector<double>& pos)
-// {
-//   size_t N = pos.size() / 2;
-
-//   std::vector<double> output(N*3);
-
-//   for ( size_t i = 0 ; i < N ; ++i ) {
-//     output[i*3+0] = pos[i*2+0];
-//     output[i*3+1] = pos[i*2+1];
-//     output[i*3+2] = std::pow(pos[i*2+0]*pos[i*2+1],.5);
-//   }
-
-//   return output;
-// }
-
-// // ----------------
-// // Python interface
-// // ----------------
-
-// namespace py = pybind11;
-
-// // wrap C++ function with NumPy array IO
-// py::array py_length(py::array_t<double, py::array::c_style | py::array::forcecast> array)
-// {
-//   // check input dimensions
-//   if ( array.ndim()     != 2 )
-//     throw std::runtime_error("Input should be 2-D NumPy array");
-//   if ( array.shape()[1] != 2 )
-//     throw std::runtime_error("Input should have size [N,2]");
-
-//   // allocate std::vector (to pass to the C++ function)
-//   std::vector<double> pos(array.size());
-
-//   // copy py::array -> std::vector
-//   std::memcpy(pos.data(),array.data(),array.size()*sizeof(double));
-
-//   // call pure C++ function
-//   std::vector<double> result = length(pos);
-
-//   ssize_t              ndim    = 2;
-//   std::vector<ssize_t> shape   = { array.shape()[0] , 3 };
-//   std::vector<ssize_t> strides = { sizeof(double)*3 , sizeof(double) };
-
-//   // return 2-D NumPy array
-//   return py::array(py::buffer_info(
-//     result.data(),                           /* data as contiguous array  */
-//     sizeof(double),                          /* size of one scalar        */
-//     py::format_descriptor<double>::format(), /* data type                 */
-//     ndim,                                    /* number of dimensions      */
-//     shape,                                   /* shape of the matrix       */
-//     strides                                  /* strides for each axis     */
-//   ));
-// }
-
-// // wrap as Python module
-// PYBIND11_MODULE(example,m)
-// {
-//   m.doc() = "pybind11 example plugin";
-
-//   m.def("length", &py_length, "Calculate the length of an array of vectors");
-// }
-
-
-
-
-  // int ScattCoeffs(const unsigned int L, const int pl, std::vector<double>& x, std::vector<std::complex<double> >& m,
-  //                 const int nmax, std::vector<std::complex<double> >& an, std::vector<std::complex<double> >& bn) {
-
-    
-// std::tuple<
-//   int, // used terms
-//   py::array_t< std::complex<double>, py::array::c_style | py::array::forcecast>,  //an
-//   py::array_t< std::complex<double>, py::array::c_style | py::array::forcecast>   //bn
-//   >
-
-py::tuple py_ScattCoeffs(
-               py::array_t<double, py::array::c_style | py::array::forcecast> py_x,
-               py::array_t< std::complex<double>, py::array::c_style | py::array::forcecast> py_m// ,
-               // const int nmax=-1, const int pl=-1
-               ) {
-  
-  std::vector<double> c_x(py_x.size());
-  std::vector< std::complex<double> > c_m(py_m.size());
-  int L = py_x.size();
-  int nmax = -1, pl = -1;
-  
-  std::memcpy(c_x.data(), py_x.data(), py_x.size()*sizeof(double));
-  std::memcpy(c_m.data(), py_m.data(), py_m.size()*sizeof( std::complex<double>));
-
-  int terms = 0;
-  std::vector<std::complex<double> > c_an, c_bn;
-  terms = nmie::ScattCoeffs( L, pl, c_x, c_m, nmax, c_an, c_bn);
-
-  auto py_an = py::array_t< std::complex<double>>(c_an.size());
-  auto py_bn = py::array_t< std::complex<double>>(c_bn.size());
-  auto py_an_buffer = py_an.request();
-  auto py_bn_buffer = py_bn.request();
-  std::complex<double> *py_an_ptr    = ( std::complex<double> *) py_an_buffer.ptr;
-  std::complex<double> *py_bn_ptr    = ( std::complex<double> *) py_bn_buffer.ptr;
-
-  std::memcpy(py_an_ptr,c_an.data(),c_an.size()*sizeof( std::complex<double>));
-  std::memcpy(py_bn_ptr,c_bn.data(),c_bn.size()*sizeof( std::complex<double>));
-
-  return py::make_tuple(terms, py_an, py_bn);
-}
-
-PYBIND11_MODULE(example, m) {
-    m.doc() = "pybind11 example plugin"; // optional module docstring
-
-    m.def("scattcoeffs", &py_ScattCoeffs, "test");
-    m.def("scattcoeffs", py::vectorize(&py_ScattCoeffs), "test");
-}