Bläddra i källkod

prepare for the switch in Python from functional C-like calls to C++ OOP

Konstantin Ladutenko 3 år sedan
förälder
incheckning
244b0cf621

+ 2 - 1
examples/Opt.Commun.-2010-Geints-Nanojets_of_dielectric_microspheres/fig1.py

@@ -40,6 +40,7 @@ import inspect
 print("Using scattnlay from ", inspect.getfile(scattnlay))
 
 npts = 351
+# npts = 51
 factor = 3 # plot extent compared to sphere radius
 index_H2O = 1.33+(1e-6)*1j
 
@@ -62,7 +63,7 @@ print("x =", x)
 print("m =", m)
 terms, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2 = scattnlay(
     np.array([x]), np.array([m]))
-print("Qsca = " + str(Qsca)+" terms = "+str(terms))
+print("   Qsca = " + str(Qsca)+" terms = "+str(terms))
 terms, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2 = scattnlay(
     np.array([x]), np.array([m]), mp=True)
 print("mp Qsca = " + str(Qsca)+" terms = "+str(terms))

+ 15 - 6
scattnlay/main.py

@@ -164,6 +164,20 @@ def expancoeffs(x, m, nmax=-1, pl=-1, mp=False):
     return terms, an, bn, cn, dn
 #expancoeffs()
 
+def scattnlay_(x, m, theta=np.zeros(0, dtype=float), nmax=-1, pl=-1, mp=False):
+    if mp:
+        from scattnlay_mp import scattnlay as scattnlay2_
+        return scattnlay2_(x, m, theta, nmax=nmax, pl=pl)
+
+    from scattnlay_dp import scattnlay as scattnlay2_
+    return scattnlay2_(x, m, theta, nmax=nmax, pl=pl)
+        # from scattnlay_dp import mie_dp as mie_
+    # mie = mie_()
+    # a = mie.GetPECLayer()
+
+    # terms[i], Qext[i], Qsca[i], Qabs[i], Qbk[i], Qpr[i], g[i], Albedo[i], S1[i], S2[i] = scattnlay_(xi, m[i], theta, nmax=nmax, pl=pl)
+    #
+    # return terms, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2
 
 def scattnlay(x, m, theta=np.zeros(0, dtype=float), nmax=-1, pl=-1, mp=False):
     """
@@ -193,11 +207,6 @@ def scattnlay(x, m, theta=np.zeros(0, dtype=float), nmax=-1, pl=-1, mp=False):
         S1, S2: Complex scattering amplitudes
     """
 
-    if mp:
-        from scattnlay_mp import scattnlay as scattnlay_
-    else:
-        from scattnlay_dp import scattnlay as scattnlay_
-
     if len(m.shape) != 1 and len(m.shape) != 2:
         raise ValueError('The relative refractive index (m) should be a 1-D or 2-D NumPy array.')
     if len(x.shape) == 1:
@@ -226,7 +235,7 @@ def scattnlay(x, m, theta=np.zeros(0, dtype=float), nmax=-1, pl=-1, mp=False):
     S2 = np.zeros((x.shape[0], theta.shape[0]), dtype=complex)
 
     for i, xi in enumerate(x):
-        terms[i], Qext[i], Qsca[i], Qabs[i], Qbk[i], Qpr[i], g[i], Albedo[i], S1[i], S2[i] = scattnlay_(xi, m[i], theta, nmax=nmax, pl=pl)
+        terms[i], Qext[i], Qsca[i], Qabs[i], Qbk[i], Qpr[i], g[i], Albedo[i], S1[i], S2[i] = scattnlay_(xi, m[i], theta, nmax=nmax, pl=pl, mp=mp)
 
     return terms, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2
 #scattnlay()

+ 1 - 1
src/nmie.cc

@@ -415,7 +415,7 @@ namespace nmie {
       ml_mie.SetFieldCoords({ConvertVector<FloatType>(Xp_vec),
 	    ConvertVector<FloatType>(Yp_vec),
 	    ConvertVector<FloatType>(Zp_vec) });
-      if (nmax != -1) ml_mie.SetMaxTerms(nmax);
+      ml_mie.SetMaxTerms(nmax);
       ml_mie.SetModeNmaxAndType(mode_n, mode_type);
 
       ml_mie.RunFieldCalculation();

+ 17 - 0
src/pb11_wrapper.cc

@@ -1,16 +1,33 @@
+#include <string>
 #include <pybind11/pybind11.h>
 #include "nmie-pybind11.hpp"
+#include "nmie.hpp"
 
 namespace py = pybind11;
 
+template<typename T>
+void declare_nmie(py::module &m, std::string &typestr) {
+  using mie_typed = nmie::MultiLayerMie<nmie::FloatType>;
+  std::string pyclass_name = std::string("mie") + typestr;
+  py::class_<mie_typed>(m, pyclass_name.c_str(), py::buffer_protocol(), py::dynamic_attr())
+      .def(py::init<>())
+      .def("GetPECLayer", &mie_typed::GetPECLayer)
+//      .def("width",     &mie_typed::width)
+//      .def("height",    &mie_typed::height)
+      ;
+}
+
 // 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<nmie::FloatType>(m, precision_name);
 
   m.def("scattcoeffs", &py_ScattCoeffs,
         "Calculate the scattering coefficients, required to calculate both the near- and far-field parameters.",