Konstantin Ladutenko 5 years ago
parent
commit
cb4529cc7a
5 changed files with 161 additions and 5 deletions
  1. 5 2
      Makefile
  2. 14 2
      src/nmie-applied-impl.hpp
  3. 2 0
      src/nmie-applied.hpp
  4. 137 0
      src/nmie-js-wrapper.cc
  5. 3 1
      src/nmie.hpp

+ 5 - 2
Makefile

@@ -42,10 +42,10 @@ standalone: scattnlay fieldnlay scattnlay-mp fieldnlay-mp
 
 # standalone programs with DP
 scattnlay: $(SRCDIR)/farfield.cc $(SRCDIR)/nmie.cc $(CXX_NMIE_HEADERS)
-	$(CXX) -DNDEBUG -O2 -Wall -std=c++11 $(SRCDIR)/farfield.cc $(SRCDIR)/nmie.cc  -lm -o scattnlay $(CXXFLAGS) $(LDFLAGS)
+	$(CXX) -DNDEBUG -O2 -Wall -std=c++11 $(SRCDIR)/farfield.cc $(SRCDIR)/nmie.cc  -lm -o scattnlay-dp $(CXXFLAGS) $(LDFLAGS)
 
 fieldnlay: $(SRCDIR)/nearfield.cc $(SRCDIR)/nmie.cc $(CXX_NMIE_HEADERS)
-	$(CXX) -DNDEBUG -O2 -Wall -std=c++11 $(SRCDIR)/nearfield.cc $(SRCDIR)/nmie.cc  -lm -o fieldnlay $(CXXFLAGS) $(LDFLAGS)
+	$(CXX) -DNDEBUG -O2 -Wall -std=c++11 $(SRCDIR)/nearfield.cc $(SRCDIR)/nmie.cc  -lm -o fieldnlay-dp $(CXXFLAGS) $(LDFLAGS)
 
 # standalone programs with MP
 scattnlay-mp: $(SRCDIR)/farfield.cc $(SRCDIR)/nmie.cc $(CXX_NMIE_HEADERS)
@@ -54,6 +54,9 @@ scattnlay-mp: $(SRCDIR)/farfield.cc $(SRCDIR)/nmie.cc $(CXX_NMIE_HEADERS)
 fieldnlay-mp: $(SRCDIR)/nearfield.cc $(SRCDIR)/nmie.cc $(CXX_NMIE_HEADERS)
 	$(CXX) -DNDEBUG -DMULTI_PRECISION=$(MULTIPREC) -O2 -Wall -std=c++11 $(SRCDIR)/nearfield.cc $(SRCDIR)/nmie.cc  -lm -o fieldnlay-mp $(CXXFLAGS) $(LDFLAGS)
 
+wasm: $(SRCDIR)/nmie-js-wrapper.cc $(CXX_NMIE_HEADERS)
+	emcc --bind -lm -Wall -O2 -std=c++11 -s WASM=1 -s NO_EXIT_RUNTIME=1 -s "EXTRA_EXPORTED_RUNTIME_METHODS=['addOnPostRun']" -o nmie.js $(SRCDIR)/nmie-js-wrapper.cc
+
 clean:
 	$(PYTHON) setup.py clean
 	$(MAKE) -f $(CURDIR)/debian/rules clean

+ 14 - 2
src/nmie-applied-impl.hpp

@@ -74,12 +74,24 @@ namespace nmie {
   // ********************************************************************** //
   template <typename FloatType>
   void MultiLayerMieApplied<FloatType>::AddTargetLayer(FloatType width, std::complex<FloatType> layer_index) {
+      this->MarkUncalculated();
+      if (width <= 0)
+          throw std::invalid_argument("Layer width should be positive!");
+      target_width_.push_back(width);
+      target_index_.push_back(layer_index);
+  }  // end of void  MultiLayerMieApplied<FloatType>::AddTargetLayer(...)
+  // ********************************************************************** //
+  // ********************************************************************** //
+  // ********************************************************************** //
+  template <typename FloatType>
+  void MultiLayerMieApplied<FloatType>::AddTargetLayerReIm(FloatType width,
+                                                           FloatType re_layer_index, FloatType im_layer_index) {
     this->MarkUncalculated();
     if (width <= 0)
       throw std::invalid_argument("Layer width should be positive!");
     target_width_.push_back(width);
-    target_index_.push_back(layer_index);
-  }  // end of void  MultiLayerMieApplied<FloatType>::AddTargetLayer(...)  
+    target_index_.push_back(std::complex<FloatType>(re_layer_index,im_layer_index));
+  }  // end of void  MultiLayerMieApplied<FloatType>::AddTargetLayer(...)
   // ********************************************************************** //
   // ********************************************************************** //
   // ********************************************************************** //

+ 2 - 0
src/nmie-applied.hpp

@@ -70,6 +70,8 @@ namespace nmie {
     // For many runs it can be convenient to separate target and coating layers.
     // Per layer
     void AddTargetLayer(FloatType layer_width, std::complex<FloatType> layer_index);
+    void AddTargetLayerReIm(FloatType layer_width,
+            FloatType re_layer_index, FloatType im_layer_index);
     void AddCoatingLayer(FloatType layer_width, std::complex<FloatType> layer_index);
     // For all layers
     void SetTargetWidth(std::vector<FloatType> width);

+ 137 - 0
src/nmie-js-wrapper.cc

@@ -0,0 +1,137 @@
+//**********************************************************************************//
+//    Copyright (C) 2009-2019  Ovidio Pena <ovidio@bytesfall.com>                   //
+//    Copyright (C) 2013-2019  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/>.         //
+//                                                                                  //
+//    @brief  Wrapper to JS
+//                                                                                  //
+//**********************************************************************************//
+#include "nmie-applied.hpp"
+#include "nmie-applied-impl.hpp"
+#include "nmie-precision.hpp"
+#include <emscripten/bind.h>
+
+using namespace emscripten;
+
+nmie::MultiLayerMieApplied<double> ml_mie;
+
+EMSCRIPTEN_BINDINGS (c) {
+        class_<nmie::MultiLayerMie<double>>("nmie_base")
+                .constructor<>()
+                .function("GetQsca",&nmie::MultiLayerMie<double>::GetQsca)
+                .function("GetQext",&nmie::MultiLayerMie<double>::GetQext)
+                .function("GetQabs",&nmie::MultiLayerMie<double>::GetQabs)
+        ;
+        class_<nmie::MultiLayerMieApplied<double>>("nmie")
+                .constructor<>()
+                .function("SetWavelength", &nmie::MultiLayerMieApplied<double>::SetWavelength)
+                .function("AddTargetLayerReIm",&nmie::MultiLayerMieApplied<double>::AddTargetLayerReIm)
+                .function("RunMieCalculation",&nmie::MultiLayerMieApplied<double>::RunMieCalculation)
+//                .function("bf",&nmie::MultiLayerMieApplied<double>::bf)
+        ;
+}
+
+//namespace nmie {
+  //**********************************************************************************//
+  // This function emulates a C call to calculate the actual scattering parameters    //
+  // and amplitudes.                                                                  //
+  //                                                                                  //
+  // Input parameters:                                                                //
+  //   L: Number of layers                                                            //
+  //   pl: Index of PEC layer. If there is none just send -1                          //
+  //   x: Array containing the size parameters of the layers [0..L-1]                 //
+  //   m: Array containing the relative refractive indexes of the layers [0..L-1]     //
+  //   nTheta: Number of scattering angles                                            //
+  //   Theta: Array containing all the scattering angles where the scattering         //
+  //          amplitudes will be calculated                                           //
+  //   nmax: Maximum number of multipolar expansion terms to be used for the          //
+  //         calculations. Only use it if you know what you are doing, otherwise      //
+  //         set this parameter to -1 and the function will calculate it              //
+  //                                                                                  //
+  // Output parameters:                                                               //
+  //   Qext: Efficiency factor for extinction                                         //
+  //   Qsca: Efficiency factor for scattering                                         //
+  //   Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca)                    //
+  //   Qbk: Efficiency factor for backscattering                                      //
+  //   Qpr: Efficiency factor for the radiation pressure                              //
+  //   g: Asymmetry factor (g = (Qext-Qpr)/Qsca)                                      //
+  //   Albedo: Single scattering albedo (Albedo = Qsca/Qext)                          //
+  //   S1, S2: Complex scattering amplitudes                                          //
+  //                                                                                  //
+  // Return value:                                                                    //
+  //   Number of multipolar expansion terms used for the calculations                 //
+  //**********************************************************************************//
+//  int nMieApplied(const unsigned int L, const int pl, std::vector<double>& x, std::vector<std::complex<double> >& m, const unsigned int nTheta, std::vector<double>& Theta, const int nmax, double *Qext, double *Qsca, double *Qabs, double *Qbk, double *Qpr, double *g, double *Albedo, std::vector<std::complex<double> >& S1, std::vector<std::complex<double> >& S2) {
+//
+//    if (x.size() != L || m.size() != L)
+//        throw std::invalid_argument("Declared number of layers do not fit x and m!");
+//    if (Theta.size() != nTheta)
+//        throw std::invalid_argument("Declared number of sample for Theta is not correct!");
+//    try {
+//      MultiLayerMieApplied<FloatType> ml_mie;
+//      ml_mie.SetLayersSize(ConvertVector<FloatType>(x));
+//      ml_mie.SetLayersIndex(ConvertComplexVector<FloatType>(m));
+//      ml_mie.SetAngles(ConvertVector<FloatType>(Theta));
+//      ml_mie.SetPECLayer(pl);
+//      ml_mie.SetMaxTerms(nmax);
+//
+//      ml_mie.RunMieCalculation();
+//
+//      *Qext = static_cast<double>(ml_mie.GetQext());
+//      *Qsca = static_cast<double>(ml_mie.GetQsca());
+//      *Qabs = static_cast<double>(ml_mie.GetQabs());
+//      *Qbk = static_cast<double>(ml_mie.GetQbk());
+//      *Qpr = static_cast<double>(ml_mie.GetQpr());
+//      *g = static_cast<double>(ml_mie.GetAsymmetryFactor());
+//      *Albedo = static_cast<double>(ml_mie.GetAlbedo());
+//      S1 = ConvertComplexVector<double>(ml_mie.GetS1());
+//      S2 = ConvertComplexVector<double>(ml_mie.GetS2());
+//
+//      return ml_mie.GetMaxTerms();
+//    } catch(const std::invalid_argument& ia) {
+//      // Will catch if  ml_mie fails or other errors.
+//      std::cerr << "Invalid argument: " << ia.what() << std::endl;
+//      throw std::invalid_argument(ia);
+//      return -1;
+//    }
+//    return 0;
+//  }
+//  int nMieApplied(const unsigned int L, std::vector<double>& x, std::vector<std::complex<double> >& m, const unsigned int nTheta, std::vector<double>& Theta, double *Qext, double *Qsca, double *Qabs, double *Qbk, double *Qpr, double *g, double *Albedo, std::vector<std::complex<double> >& S1, std::vector<std::complex<double> >& S2) {
+//    return nmie::nMieApplied(L, -1, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
+//  }
+//  int nMieApplied(const unsigned int L, const int pl, std::vector<double>& x, std::vector<std::complex<double> >& m, const unsigned int nTheta, std::vector<double>& Theta, double *Qext, double *Qsca, double *Qabs, double *Qbk, double *Qpr, double *g, double *Albedo, std::vector<std::complex<double> >& S1, std::vector<std::complex<double> >& S2) {
+//    return nmie::nMieApplied(L, pl, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
+//  }
+//  int nMieApplied(const unsigned int L, std::vector<double>& x, std::vector<std::complex<double> >& m, const unsigned int nTheta, std::vector<double>& Theta, const int nmax, double *Qext, double *Qsca, double *Qabs, double *Qbk, double *Qpr, double *g, double *Albedo, std::vector<std::complex<double> >& S1, std::vector<std::complex<double> >& S2) {
+//    return nmie::nMieApplied(L, -1, x, m, nTheta, Theta, nmax, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
+//  }
+  // ********************************************************************** //
+  // ********************************************************************** //
+  // ********************************************************************** //
+
+//}  // end of namespace nmie

+ 3 - 1
src/nmie.hpp

@@ -37,7 +37,9 @@
 #include <cstdlib>
 #include <iostream>
 #include <vector>
+#ifdef MULTI_PRECISION
 #include <boost/math/constants/constants.hpp>
+#endif
 namespace nmie {
   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);
   int nMie(const unsigned int L, const int pl, std::vector<double>& x, std::vector<std::complex<double> >& m, const unsigned int nTheta, std::vector<double>& Theta, const int nmax, double *Qext, double *Qsca, double *Qabs, double *Qbk, double *Qpr, double *g, double *Albedo, std::vector<std::complex<double> >& S1, std::vector<std::complex<double> >& S2);
@@ -88,7 +90,7 @@ namespace nmie {
 
     // Set a fixed value for the maximun number of terms
     void SetMaxTerms(int nmax);
-    // Get maximun number of terms
+    // Get maximum number of terms
     int GetMaxTerms() {return nmax_;};
 
     bool isMieCalculated(){return isMieCalculated_;};