Lemma is an Electromagnetics API
Du kannst nicht mehr als 25 Themen auswählen Themen müssen mit entweder einem Buchstaben oder einer Ziffer beginnen. Sie können Bindestriche („-“) enthalten und bis zu 35 Zeichen lang sein.

pyFDEM1D.cpp 20KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341
  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");
  26. // lifecycle
  27. WireAntenna.def(py::init(&Lemma::WireAntenna::NewSP))
  28. .def_static("DeSerialize", py::overload_cast<const std::string&>(&Lemma::WireAntenna::DeSerialize),
  29. "Construct object from yaml representation")
  30. // print
  31. .def("Serialize", &Lemma::WireAntenna::Print, "YAML representation of the class")
  32. .def("__repr__", &Lemma::WireAntenna::Print)
  33. // modifiers
  34. .def("SetNumberOfPoints", &Lemma::WireAntenna::SetNumberOfPoints,
  35. "Sets the number of points comprising the antenna")
  36. .def("SetPoint", py::overload_cast<const int&, const Lemma::Real&, const Lemma::Real&, const Lemma::Real&>(&Lemma::WireAntenna::SetPoint),
  37. "Sets a point in the antenna")
  38. .def("SetPoint", py::overload_cast<const int&, const Lemma::Vector3r&>(&Lemma::WireAntenna::SetPoint),
  39. "Sets a point in the antenna")
  40. .def("SetNumberOfTurns", &Lemma::WireAntenna::SetNumberOfTurns, "Sets the number of turns of the antenna")
  41. .def("SetNumberOfFrequencies", &Lemma::WireAntenna::SetNumberOfFrequencies,
  42. "Sets the number of frequencies of the transmitter")
  43. .def("SetFrequency", &Lemma::WireAntenna::SetFrequency, "Sets a single frequency of the transmitter")
  44. .def("SetCurrent", &Lemma::WireAntenna::SetCurrent, "Sets the current of the transmitter in amps")
  45. // accessors
  46. .def("GetCurrent", &Lemma::WireAntenna::GetCurrent, "Returns the current of the transmitter in amps")
  47. .def("GetPoints", &Lemma::WireAntenna::GetPoints, "Returns the points defining the transmitter")
  48. .def("GetNumberOfDipoles", &Lemma::WireAntenna::GetNumberOfDipoles,
  49. "Returns the number of dipoles defining the transmitter")
  50. .def("GetNumberOfFrequencies", &Lemma::WireAntenna::GetNumberOfFrequencies,
  51. "Returns the number of frequencies for the transmitter")
  52. .def("IsHorizontallyPlanar", &Lemma::WireAntenna::IsHorizontallyPlanar, "Returns true if the transmitter is flat")
  53. .def("GetName", &Lemma::WireAntenna::GetName, "Returns the class name of the object")
  54. // operations
  55. .def("ApproximateWithElectricDipoles", &Lemma::WireAntenna::ApproximateWithElectricDipoles,
  56. "Approximates loop with electric dipoles")
  57. ;
  58. py::class_<Lemma::PolygonalWireAntenna, std::shared_ptr<Lemma::PolygonalWireAntenna> > PolygonalWireAntenna(m,
  59. "PolygonalWireAntenna", WireAntenna);
  60. // lifecycle
  61. PolygonalWireAntenna.def(py::init(&Lemma::PolygonalWireAntenna::NewSP))
  62. .def_static("DeSerialize", py::overload_cast<const std::string&>(&Lemma::PolygonalWireAntenna::DeSerialize),
  63. "Construct object from yaml representation")
  64. // print
  65. .def("__repr__", &Lemma::PolygonalWireAntenna::Print)
  66. .def("Serialize", &Lemma::PolygonalWireAntenna::Print, "YAML representation of the class")
  67. // accessors
  68. .def("GetName", &Lemma::PolygonalWireAntenna::GetName, "Returns the name of the class")
  69. // operations
  70. .def("ApproximateWithElectricDipoles", &Lemma::PolygonalWireAntenna::ApproximateWithElectricDipoles,
  71. "Approximates loop with series of electric dipoles around loop")
  72. // modifiers
  73. .def("SetMinDipoleRatio", &Lemma::PolygonalWireAntenna::SetMinDipoleRatio,
  74. "Sets the minimum dipole ratio use, smaller values increase precision")
  75. .def("SetMinDipoleMoment", &Lemma::PolygonalWireAntenna::SetMinDipoleMoment,
  76. "Sets the minimum dipole moment which will be used, smaller values increase precision and computational time")
  77. .def("SetMaxDipoleMoment", &Lemma::PolygonalWireAntenna::SetMaxDipoleMoment,
  78. "Sets the maximum dipole moment which will be used, smaller values increase precision and computational time")
  79. ;
  80. py::class_<Lemma::DipoleSource, std::shared_ptr<Lemma::DipoleSource> > DipoleSource(m, "DipoleSource");
  81. // lifecycle
  82. DipoleSource.def(py::init(&Lemma::DipoleSource::NewSP))
  83. .def_static("DeSerialize", py::overload_cast<const std::string&>(&Lemma::DipoleSource::DeSerialize),
  84. "Construct object from yaml representation")
  85. // print
  86. .def("Serialize", &Lemma::DipoleSource::Print, "YAML representation of the class")
  87. .def("__repr__", &Lemma::DipoleSource::Print)
  88. // accessors
  89. .def("GetName", &Lemma::DipoleSource::GetName, "Returns the name of the class")
  90. .def("GetNumberOfFrequencies", &Lemma::DipoleSource::GetNumberOfFrequencies,
  91. "Returns the number of frequencies")
  92. .def("GetFrequencies", &Lemma::DipoleSource::GetFrequencies, "Returns an array of frequencies")
  93. .def("GetFrequency", &Lemma::DipoleSource::GetFrequency, "Returns the frequency of the argument index")
  94. .def("GetAngularFrequency", &Lemma::DipoleSource::GetAngularFrequency,
  95. "Returns the angular frequency of the argument index")
  96. .def("GetPhase", &Lemma::DipoleSource::GetPhase, "Returns the phase of the dipole")
  97. .def("GetMoment", &Lemma::DipoleSource::GetMoment, "Returns the moment of the dipole")
  98. .def("GetLocation", py::overload_cast< >(&Lemma::DipoleSource::GetLocation), "Returns the location of the dipole")
  99. .def("GetPolarisation", &Lemma::DipoleSource::GetPolarisation, "Returns the polarisation of the dipole")
  100. // modifiers
  101. .def("SetLocation", py::overload_cast<const Lemma::Vector3r&> (&Lemma::DipoleSource::SetLocation),
  102. "Sets the location of the dipole")
  103. .def("SetPolarisation", py::overload_cast<const Lemma::Vector3r&> (&Lemma::DipoleSource::SetPolarisation),
  104. "Sets the polarisation of the dipole")
  105. .def("SetType", &Lemma::DipoleSource::SetType, "Sets the type")
  106. .def("SetMoment", &Lemma::DipoleSource::SetMoment, "Sets the moment of the dipole")
  107. .def("SetPhase", &Lemma::DipoleSource::SetPhase, "Sets the phase of the dipole")
  108. .def("SetNumberOfFrequencies", &Lemma::DipoleSource::SetNumberOfFrequencies,
  109. "Sets the number of frequencies to calculate for the dipole")
  110. .def("SetFrequency", &Lemma::DipoleSource::SetFrequency,
  111. "Sets a single frequency, first argument is index, second argument is frequency")
  112. .def("SetFrequencies", &Lemma::DipoleSource::SetFrequencies,
  113. "Sets all frequencies, argument is numpy array of frequencies")
  114. ;
  115. py::class_<Lemma::LayeredEarthEM, std::shared_ptr<Lemma::LayeredEarthEM>, Lemma::EarthModel > LayeredEarthEM(m, "LayeredEarthEM");
  116. // lifecycle
  117. LayeredEarthEM.def(py::init(&Lemma::LayeredEarthEM::NewSP))
  118. .def_static("DeSerialize", py::overload_cast<const std::string&>
  119. (&Lemma::LayeredEarthEM::DeSerialize),"Construct object from yaml representation")
  120. // print
  121. .def("Serialize", &Lemma::LayeredEarthEM::Print, "YAML representation of the class")
  122. .def("__repr__", &Lemma::LayeredEarthEM::Print)
  123. // accessors
  124. .def("GetName", &Lemma::LayeredEarthEM::GetName, "Returns the name of the class")
  125. .def("GetLayerConductivity", py::overload_cast<>(&Lemma::LayeredEarthEM::GetLayerConductivity),
  126. "Returns the conductivity of all layers")
  127. .def("GetLayerConductivity1", py::overload_cast<const int&>(&Lemma::LayeredEarthEM::GetLayerConductivity),
  128. "Returns the conductivity of the specified layer")
  129. .def("GetLayerSusceptibility", py::overload_cast<>(&Lemma::LayeredEarthEM::GetLayerSusceptibility),
  130. "Returns the susceptibility of all layers")
  131. .def("GetLayerSusceptibility1", py::overload_cast<const int&>(&Lemma::LayeredEarthEM::GetLayerSusceptibility),
  132. "Returns the susceptibilty of the specified layer")
  133. .def("GetLayerLowFreqSusceptibility", py::overload_cast<>(&Lemma::LayeredEarthEM::GetLayerLowFreqSusceptibility),
  134. "Returns the low frequqncy permitivity of all layers")
  135. .def("GetLayerLowFreqSusceptibility1", py::overload_cast<const int&>(&Lemma::LayeredEarthEM::GetLayerLowFreqSusceptibility),
  136. "Returns the low frequency permitivity of the specified layer")
  137. .def("GetLayerHighFreqSusceptibility", py::overload_cast<>(&Lemma::LayeredEarthEM::GetLayerHighFreqSusceptibility),
  138. "Returns the low frequency permitivity of all layers")
  139. .def("GetLayerHighFreqSusceptibility1", py::overload_cast<const int&>(&Lemma::LayeredEarthEM::GetLayerHighFreqSusceptibility),
  140. "Returns the low frequency permitivity of the specified layer")
  141. .def("GetLayerTauSusceptibility", py::overload_cast<>(&Lemma::LayeredEarthEM::GetLayerTauSusceptibility),
  142. "Returns the tau permitivity of all layers")
  143. .def("GetLayerTauSusceptibility1", py::overload_cast<const int&>(&Lemma::LayeredEarthEM::GetLayerTauSusceptibility),
  144. "Returns the tau permitivity of the specified layer")
  145. .def("GetLayerBreathSusceptibility", py::overload_cast<>(&Lemma::LayeredEarthEM::GetLayerBreathSusceptibility),
  146. "Returns the breth permitivity of all layers")
  147. .def("GetLayerBreathSusceptibility1", py::overload_cast<const int&>(&Lemma::LayeredEarthEM::GetLayerBreathSusceptibility),
  148. "Returns the breath permitivity of the specified layer")
  149. .def("GetLayerPermitivity", py::overload_cast<>(&Lemma::LayeredEarthEM::GetLayerPermitivity),
  150. "Returns the permitivity of all layers")
  151. .def("GetLayerPermitivity1", py::overload_cast<const int&>(&Lemma::LayeredEarthEM::GetLayerPermitivity),
  152. "Returns the permitivity of the specified layer")
  153. .def("GetLayerLowFreqPermitivity", py::overload_cast<>(&Lemma::LayeredEarthEM::GetLayerLowFreqPermitivity),
  154. "Returns the low frequqncy permitivity of all layers")
  155. .def("GetLayerLowFreqPermitivity1", py::overload_cast<const int&>(&Lemma::LayeredEarthEM::GetLayerLowFreqPermitivity),
  156. "Returns the low frequency permitivity of the specified layer")
  157. .def("GetLayerHighFreqPermitivity", py::overload_cast<>(&Lemma::LayeredEarthEM::GetLayerHighFreqPermitivity),
  158. "Returns the low frequency permitivity of all layers")
  159. .def("GetLayerHighFreqPermitivity1", py::overload_cast<const int&>(&Lemma::LayeredEarthEM::GetLayerHighFreqPermitivity),
  160. "Returns the low frequency permitivity of the specified layer")
  161. .def("GetLayerTauPermitivity", py::overload_cast<>(&Lemma::LayeredEarthEM::GetLayerTauPermitivity),
  162. "Returns the tau permitivity of all layers")
  163. .def("GetLayerTauPermitivity1", py::overload_cast<const int&>(&Lemma::LayeredEarthEM::GetLayerTauPermitivity),
  164. "Returns the tau permitivity of the specified layer")
  165. .def("GetLayerBreathPermitivity", py::overload_cast<>(&Lemma::LayeredEarthEM::GetLayerBreathPermitivity),
  166. "Returns the breth permitivity of all layers")
  167. .def("GetLayerBreathPermitivity1", py::overload_cast<const int&>(&Lemma::LayeredEarthEM::GetLayerBreathPermitivity),
  168. "Returns the breath permitivity of the specified layer")
  169. // modifiers
  170. .def("SetNumberOfLayers", &Lemma::LayeredEarthEM::SetNumberOfLayers,
  171. "Sets the number of layers in the model")
  172. .def("SetLayerConductivity", py::overload_cast< const Lemma::VectorXcr& >(&Lemma::LayeredEarthEM::SetLayerConductivity),
  173. "Sets the conductivity of the layers, the input is a complex array of conductivity")
  174. .def("SetLayerConductivity1", py::overload_cast< const int&, const Lemma::Complex& >(&Lemma::LayeredEarthEM::SetLayerConductivity),
  175. "Sets the conductivity of a single layer, the first input is the layer index, and the secondinput is a complex conductivity")
  176. .def("SetLayerThickness", &Lemma::LayeredEarthEM::SetLayerThickness,
  177. "Sets the thickness of layers, excluding the air and bottom which are infinite")
  178. .def("SetLayerHighFreqSusceptibility", &Lemma::LayeredEarthEM::SetLayerHighFreqSusceptibility,
  179. "Sets the high frequency susceptibility for Cole-COle model")
  180. .def("SetLayerLowFreqSusceptibility", &Lemma::LayeredEarthEM::SetLayerLowFreqSusceptibility,
  181. "Sets the low frequency susceptibility for Cole-COle model")
  182. .def("SetLayerBreathSusceptibility", &Lemma::LayeredEarthEM::SetLayerBreathSusceptibility,
  183. "Sets thesusceptibility breath for Cole-COle model")
  184. .def("SetLayerHighFreqPermitivity", &Lemma::LayeredEarthEM::SetLayerHighFreqPermitivity,
  185. "Sets the high frequency permitivity for Cole-COle model")
  186. .def("SetLayerLowFreqPermitivity", &Lemma::LayeredEarthEM::SetLayerLowFreqPermitivity,
  187. "Sets the low frequency permitivity for Cole-COle model")
  188. .def("SetLayerBreathPermitivity", &Lemma::LayeredEarthEM::SetLayerBreathPermitivity,
  189. "Sets the permitivity breath for Cole-COle model")
  190. // methods
  191. .def("EvaluateColeColeModel", &Lemma::LayeredEarthEM::EvaluateColeColeModel,
  192. "Calculates complex resistivity based on cole-cole parameters")
  193. ;
  194. py::class_<Lemma::EMEarth1D, std::shared_ptr<Lemma::EMEarth1D> >
  195. EMEarth1D(m, "EMEarth1D");
  196. // lifecycle
  197. EMEarth1D.def(py::init(&Lemma::EMEarth1D::NewSP))
  198. //.def_static("DeSerialize", py::overload_cast<const std::string&>
  199. // (&Lemma::EMEarth1D::DeSerialize),"Construct object from yaml representation")
  200. // print
  201. .def("Serialize", &Lemma::EMEarth1D::Print, "YAML representation of the class")
  202. .def("__repr__", &Lemma::EMEarth1D::Print)
  203. // accessors
  204. .def("GetName", &Lemma::EMEarth1D::GetName, "Returns the name of the class")
  205. .def("GetFieldPoints", &Lemma::EMEarth1D::GetFieldPoints, "Returns the FieldPoint class")
  206. // modifiers
  207. .def("AttachWireAntenna", &Lemma::EMEarth1D::AttachWireAntenna,
  208. "Sets the wire antenna to use for calculations")
  209. .def("AttachDipoleSource", &Lemma::EMEarth1D::AttachDipoleSource,
  210. "Sets a DipoleSource to use for calculations")
  211. .def("AttachFieldPoints", &Lemma::EMEarth1D::AttachFieldPoints,
  212. "Sets the FieldPoints to use for calculations")
  213. .def("AttachLayeredEarthEM", &Lemma::EMEarth1D::AttachLayeredEarthEM,
  214. "Sets the LayeredEarthEM to use for calculations")
  215. .def("SetFieldsToCalculate", &Lemma::EMEarth1D::SetFieldsToCalculate,
  216. "Sets which fields to calculate")
  217. .def("SetHankelTransformMethod", &Lemma::EMEarth1D::SetHankelTransformMethod,
  218. "Sets which Hankel transform to use")
  219. .def("SetTxRxMode", &Lemma::EMEarth1D::SetTxRxMode,
  220. "Sets the TxRx mode flag")
  221. //methods
  222. #ifdef KIHALEE_EM1D
  223. .def("MakeCalc", &Lemma::EMEarth1D::MakeCalc, "Calls KiHa Lee's EM1D FORTRAN77 code")
  224. #endif
  225. .def("MakeCalc3", &Lemma::EMEarth1D::MakeCalc3, "Native Lemma EM calculation")
  226. .def("CalculateWireAntennaFields", &Lemma::EMEarth1D::CalculateWireAntennaFields,
  227. "Native Lemma calculation of a wire antenna")
  228. ;
  229. py::class_<Lemma::FieldPoints, std::shared_ptr<Lemma::FieldPoints> >
  230. FieldPoints(m, "FieldPoints");
  231. // lifecycle
  232. FieldPoints.def(py::init(&Lemma::FieldPoints::NewSP))
  233. .def_static("DeSerialize", py::overload_cast<const std::string&>
  234. (&Lemma::FieldPoints::DeSerialize),"Construct object from yaml representation")
  235. // print
  236. .def("Serialize", &Lemma::FieldPoints::Print, "YAML representation of the class")
  237. .def("__repr__", &Lemma::FieldPoints::Print)
  238. // modifiers
  239. .def("SetNumberOfPoints", &Lemma::FieldPoints::SetNumberOfPoints,
  240. "Sets the number of locations to make calculations on.")
  241. .def("SetLocation", py::overload_cast< const int&, const Lemma::Vector3r& >
  242. (&Lemma::FieldPoints::SetLocation), "Sets the location of the index-specified point." )
  243. .def("SetLocation", py::overload_cast< const int&,
  244. const Lemma::Real&, const Lemma::Real&, const Lemma::Real& >
  245. (&Lemma::FieldPoints::SetLocation),
  246. "Sets the location of the index-specified point with the three coordinates")
  247. // accessors
  248. .def("GetNumberOfPoints", &Lemma::FieldPoints::GetNumberOfPoints,
  249. "Returns the number of locations to make calculations on.")
  250. .def("GetLocations", &Lemma::FieldPoints::GetLocations,
  251. "Returns the locations which calculations are made on.")
  252. .def("GetLocationsMat", &Lemma::FieldPoints::GetLocationsMat,
  253. "Returns a matrix of the locations which calculations are made on.")
  254. .def("GetLocation", &Lemma::FieldPoints::GetLocation,
  255. "Returns the location of the specified index.")
  256. .def("GetLocationX", &Lemma::FieldPoints::GetLocationX,
  257. "Returns the northing (x) location of the specified index.")
  258. .def("GetLocationY", &Lemma::FieldPoints::GetLocationY,
  259. "Returns the easting (y) location of the specified index.")
  260. .def("GetLocationZ", &Lemma::FieldPoints::GetLocationZ,
  261. "Returns the altitude/depth (z) location of the specified index.")
  262. .def("GetEfield", py::overload_cast< > (&Lemma::FieldPoints::GetEfield),
  263. "Returns the electric field for all frequencies.")
  264. .def("GetEfield", py::overload_cast< const int& > (&Lemma::FieldPoints::GetEfield),
  265. "Returns the electric field for the specified frequency index.")
  266. .def("GetEfield", py::overload_cast< const int&, const int& > (&Lemma::FieldPoints::GetEfield),
  267. "Returns the electric field for the specified frequency and location index.")
  268. .def("GetEfieldMat", &Lemma::FieldPoints::GetEfieldMat,
  269. "Returns the electric field for the specified frequency.")
  270. .def("GetHfield", py::overload_cast< > (&Lemma::FieldPoints::GetHfield),
  271. "Returns the H field for all frequencies.")
  272. .def("GetHfield", py::overload_cast< const int& > (&Lemma::FieldPoints::GetHfield),
  273. "Returns the H field for the specified frequency index.")
  274. .def("GetHfield", py::overload_cast< const int&, const int& > (&Lemma::FieldPoints::GetHfield),
  275. "Returns the H field for the specified frequency and location index.")
  276. //.def("GetBfield", py::overload_cast< const int&, const int& > (&Lemma::FieldPoints::GetBfield),
  277. // "Returns the magnetic (B) field for the specified frequency and location index.")
  278. .def("GetHfieldMat", &Lemma::FieldPoints::GetHfieldMat,
  279. "Returns the H field for the specified frequency.")
  280. .def("GetMask", &Lemma::FieldPoints::MaskPoint, "Return the mask boolean value for the specified index")
  281. // methods
  282. .def("ClearFields", &Lemma::FieldPoints::ClearFields, "Clears calculated fields")
  283. .def("MaskPoint", &Lemma::FieldPoints::MaskPoint, "Masks the index resulting in no calculation")
  284. .def("UnMaskPoint", &Lemma::FieldPoints::UnMaskPoint, "Unmasks the index resulting in a calculation")
  285. ;
  286. }