Main Lemma Repository
You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108
  1. /* This file is part of Lemma, a geophysical modelling and inversion API.
  2. * More information is available at http://lemmasoftware.org
  3. */
  4. /* This Source Code Form is subject to the terms of the Mozilla Public
  5. * License, v. 2.0. If a copy of the MPL was not distributed with this
  6. * file, You can obtain one at http://mozilla.org/MPL/2.0/.
  7. */
  8. /**
  9. * @file
  10. * @date 22/04/19 14:06:32
  11. * @version $Id$
  12. * @author Trevor Irons (ti)
  13. * @email Trevor.Irons@utah.edu
  14. * @copyright Copyright (c) 2019, University of Utah
  15. * @copyright Copyright (c) 2019, Lemma Software, LLC
  16. */
  17. #include <pybind11/pybind11.h>
  18. #include <pybind11/iostream.h>
  19. #include <pybind11/eigen.h>
  20. #include "FDEM1D"
  21. namespace py = pybind11;
  22. PYBIND11_MODULE(FDEM1D, m) {
  23. py::add_ostream_redirect(m, "ostream_redirect");
  24. m.doc() = "Python binding of Lemma::FDEM1D, additional details can be found at https://lemmasoftware.org";
  25. py::class_<Lemma::WireAntenna, std::shared_ptr<Lemma::WireAntenna> > WireAntenna(m, "WireAntenna"); //, py::dynamic_attr());
  26. // lifecycle
  27. WireAntenna.def(py::init(&Lemma::WireAntenna::NewSP))
  28. // print
  29. .def("__repr__", &Lemma::WireAntenna::Print)
  30. .def("Serialize", &Lemma::WireAntenna::Print, "YAML representation of the class")
  31. .def_static("DeSerialize", py::overload_cast<const std::string&>(&Lemma::WireAntenna::DeSerialize), "Construct object from yaml representation")
  32. // modifiers
  33. .def("SetNumberOfPoints", &Lemma::WireAntenna::SetNumberOfPoints, "Sets the number of points comprising the antenna")
  34. .def("SetPoint", py::overload_cast<const int&, const Lemma::Real&, const Lemma::Real&, const Lemma::Real&>(&Lemma::WireAntenna::SetPoint), "Sets a point in the antenna")
  35. .def("SetPoint", py::overload_cast<const int&, const Lemma::Vector3r&>(&Lemma::WireAntenna::SetPoint), "Sets a point in the antenna")
  36. .def("SetNumberOfTurns", &Lemma::WireAntenna::SetNumberOfTurns, "Sets the number of turns of the antenna")
  37. .def("SetNumberOfFrequencies", &Lemma::WireAntenna::SetNumberOfFrequencies, "Sets the number of frequencies of the transmitter")
  38. .def("SetFrequency", &Lemma::WireAntenna::SetFrequency, "Sets a single frequency of the transmitter")
  39. .def("SetCurrent", &Lemma::WireAntenna::SetCurrent, "Sets the current of the transmitter in amps")
  40. // accessors
  41. .def("GetCurrent", &Lemma::WireAntenna::GetCurrent, "Returns the current of the transmitter in amps")
  42. .def("GetPoints", &Lemma::WireAntenna::GetPoints, "Returns the points defining the transmitter")
  43. .def("GetNumberOfDipoles", &Lemma::WireAntenna::GetNumberOfDipoles, "Returns the number of dipoles defining the transmitter")
  44. .def("GetNumberOfFrequencies", &Lemma::WireAntenna::GetNumberOfFrequencies, "Returns the number of frequencies for the transmitter")
  45. .def("IsHorizontallyPlanar", &Lemma::WireAntenna::IsHorizontallyPlanar, "Returns true if the transmitter is flat")
  46. .def("GetName", &Lemma::WireAntenna::GetName, "Returns the class name of the object")
  47. // operations
  48. .def("ApproximateWithElectricDipoles", &Lemma::WireAntenna::ApproximateWithElectricDipoles, "Approximates loop with electric dipoles")
  49. ;
  50. py::class_<Lemma::PolygonalWireAntenna, std::shared_ptr<Lemma::PolygonalWireAntenna> > PolygonalWireAntenna(m, "PolygonalWireAntenna", WireAntenna); //, py::dynamic_attr())
  51. // lifecycle
  52. PolygonalWireAntenna.def(py::init(&Lemma::PolygonalWireAntenna::NewSP))
  53. // print
  54. .def("__repr__", &Lemma::PolygonalWireAntenna::Print)
  55. .def("Serialize", &Lemma::PolygonalWireAntenna::Print, "YAML representation of the class")
  56. // accessors
  57. .def("GetName", &Lemma::PolygonalWireAntenna::GetName, "Returns the name of the class")
  58. // modifiers
  59. //.def("SetNumberOfPoints", &Lemma::WireAntenna::SetNumberOfPoints, "Sets the number of points comprising the antenna")
  60. // operations
  61. .def("ApproximateWithElectricDipoles", &Lemma::PolygonalWireAntenna::ApproximateWithElectricDipoles, "Approximates loop with series of electric dipoles around loop")
  62. ;
  63. }
  64. /*
  65. BOOST_PYTHON_MODULE(pyLemmaCore)
  66. {
  67. Py_Initialize();
  68. np::initialize();
  69. bp::class_<Lemma::PolygonalWireAntenna, boost::noncopyable> ("PolygonalWireAntenna", boost::python::no_init)
  70. // print
  71. .def(boost::python::self_ns::str(bp::self))
  72. // Lifecycle
  73. .def("__init__", boost::python::make_constructor(&Lemma::PolygonalWireAntenna::NewSP))
  74. //.def("Serialize", &Lemma::PolygonalWireAntenna::Serialize)
  75. //.def("DeSerialize", &Lemma::PolygonalWireAntenna::DeSerialize)
  76. // Accessors
  77. ;
  78. }
  79. */