Lemma is an Electromagnetics API
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.

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319
  1. // ===========================================================================
  2. //
  3. // Filename: utdipolesource.cpp
  4. //
  5. // Description:
  6. //
  7. // Version: 0.0
  8. // Created: 12/02/2009 11:57:14 AM
  9. // Revision: none
  10. // Compiler: g++ (c++)
  11. //
  12. // Author: Trevor Irons (ti)
  13. //
  14. // Organisation: Colorado School of Mines (CSM)
  15. // United States Geological Survey (USGS)
  16. //
  17. // Email: tirons@mines.edu, tirons@usgs.gov
  18. //
  19. // This program is free software: you can redistribute it and/or modify
  20. // it under the terms of the GNU General Public License as published by
  21. // the Free Software Foundation, either version 3 of the License, or
  22. // (at your option) any later version.
  23. //
  24. // This program is distributed in the hope that it will be useful,
  25. // but WITHOUT ANY WARRANTY; without even the implied warranty of
  26. // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
  27. // GNU General Public License for more details.
  28. //
  29. // You should have received a copy of the GNU General Public License
  30. // along with this program. If not, see <http://www.gnu.org/licenses/>.
  31. //
  32. // ===========================================================================
  33. #include <iostream>
  34. #include <fstream>
  35. #include "dipolesource.h"
  36. #include "layeredearth.h"
  37. #include "receiverpoints.h"
  38. #include "emearth1d.h"
  39. #include "WireAntenna.h"
  40. #include "PolygonalWireAntenna.h"
  41. #if LEMMAUSEVTK
  42. #include "vtkRenderer.h"
  43. #include "vtkRenderWindow.h"
  44. #include "vtkRenderWindowInteractor.h"
  45. #include "vtkRenderLargeImage.h"
  46. #include "vtkPNGWriter.h"
  47. #endif // ----- not LEMMA_USE_VTK -----
  48. // For testing purposes disable VTK and run scale.sh
  49. //#undef LEMMAUSEVTK
  50. #include "timer.h"
  51. using namespace Lemma;
  52. double randDouble(double low, double high) {
  53. //srand(time(0));
  54. double temp;
  55. /* swap low & high around if the user makes no sense */
  56. if (low > high) {
  57. temp = low;
  58. low = high;
  59. high = temp;
  60. }
  61. /* calculate the random number & return it */
  62. temp = (rand() / (static_cast<double>(RAND_MAX) + 1.0))
  63. * (high - low) + low;
  64. return temp;
  65. }
  66. int main() {
  67. // Keep track of time
  68. jsw_timer timer;
  69. srand(time(0));
  70. PolygonalWireAntenna *pa = PolygonalWireAntenna::New();
  71. pa->SetNumberOfFrequencies(1);
  72. pa->SetFrequency(0, 1000);
  73. pa->SetNumberOfPoints(5);
  74. pa->SetPoint(0, Vector3r( 0, 0, -.001));
  75. pa->SetPoint(1, Vector3r( 100, 0, -.001));
  76. pa->SetPoint(2, Vector3r( 100, 100, -.001));
  77. pa->SetPoint(3, Vector3r( 0, 100, -.001));
  78. pa->SetPoint(4, Vector3r( 0, 0, -.001));
  79. pa->SetCurrent(1.);
  80. pa->SetNumberOfTurns(1);
  81. //Vector3r rp = Vector3r::Random(3);
  82. //rp << 150., 10, 0.;
  83. //rp << -27.1456, 15.2350, -1e-3;
  84. //rp << randDouble(-35,35), randDouble(-35,35), randDouble(-35,35);
  85. //rp << 3.22806, -13.1548, 14.9695;
  86. //rp.setRandom(3);
  87. //std::cout << "rp " << rp.transpose() << std::endl;
  88. //pa->ApproximateWithElectricDipoles(rp);
  89. WireAntenna *wire = WireAntenna::New();
  90. wire->SetNumberOfPoints(5);
  91. wire->SetPoint(0, Vector3r( 0, 0, -1e-3));
  92. wire->SetPoint(1, Vector3r( 10, 0, -1e-3));
  93. wire->SetPoint(2, Vector3r( 10, 10, -1e-3));
  94. wire->SetPoint(3, Vector3r( 0, 10, -1e-3));
  95. wire->SetPoint(4, Vector3r( 0, 0, -1e-3));
  96. // TODO change wire antennae to use my class
  97. //wire->SetNumberOfFrequencies(1);
  98. wire->SetCurrent(1.);
  99. wire->SetNumberOfFrequencies(1);
  100. wire->SetFrequency(0, 1000);
  101. wire->SetNumberOfTurns(1);
  102. //wire->ApproximateWithElectricDipoles(5);
  103. // Define model
  104. VectorXcr sigma(2);
  105. sigma << Complex(0.,0), Complex(.1,0);
  106. VectorXr thick(1);
  107. thick << 10;
  108. LayeredEarthEM *earth = LayeredEarthEM::New();
  109. earth->SetNumberOfLayers(2);
  110. earth->SetLayerConductivity(sigma);
  111. //earth->SetLayerThickness(thick);
  112. // Receivers
  113. ReceiverPoints *receivers = ReceiverPoints::New();
  114. Vector3r loc;
  115. Real ox = 50.561 ;
  116. Real oy = 105.235 ;
  117. Real depth = -3.75e1;
  118. Real depth2 = depth;
  119. Real dx = 1.;
  120. int nz = 1;
  121. receivers->SetNumberOfReceivers(nz);
  122. int ir = 0;
  123. for (int iz=0; iz<nz; ++iz) {
  124. loc << ox, oy, depth;
  125. receivers->SetLocation(ir, loc);
  126. depth += dx;
  127. ++ ir;
  128. }
  129. // EmEarth
  130. EMEarth1D *EmEarth = EMEarth1D::New();
  131. //EmEarth->AttachWireAntenna(wire);
  132. EmEarth->AttachWireAntenna(pa);
  133. EmEarth->AttachLayeredEarthEM(earth);
  134. EmEarth->AttachReceiverPoints(receivers);
  135. EmEarth->SetFieldsToCalculate(H);
  136. //EmEarth->SetHankelTransformMethod(GAUSSIANQUADRATURE);
  137. // Do calculation
  138. timer.begin();
  139. EmEarth->CalculateWireAntennaFields();
  140. Real paTime = timer.end();
  141. std::cout << "Polygonal wire antennae time: " << paTime << "\n";
  142. //EmEarth->AttachWireAntenna(wire);
  143. //timer.begin();
  144. //EmEarth->CalculateWireAntennaFields();
  145. //Real waTime = timer.end();
  146. //std::cout << "Fixed wire antennae time: " << waTime << "\n";
  147. depth = depth2;
  148. std::fstream real("reale_lay.dat", std::ios::out);
  149. std::fstream imag("image_lay.dat", std::ios::out);
  150. std::fstream hreal("real_lay.dat", std::ios::out);
  151. std::fstream himag("imag_lay.dat", std::ios::out);
  152. for (int iz=0; iz<nz; ++iz) {
  153. Vector3cr temp = receivers->GetEfield(0,iz);
  154. real << ox << "\t" << oy << "\t" << depth << "\t"
  155. << temp(0).real() << "\t" << temp(1).real()
  156. << "\t" << temp(2).real() << std::endl;
  157. imag << ox << "\t" << oy << "\t" << depth << "\t"
  158. << std::imag(temp(0)) << "\t" << std::imag(temp(1))
  159. << "\t" << std::imag(temp(2)) << std::endl;
  160. temp = receivers->GetHfield(0, iz);
  161. hreal << ox << "\t" << oy << "\t" << depth << "\t"
  162. << std::real(temp(0)) << "\t" << std::real(temp(1))
  163. << "\t" << std::real(temp(2)) << std::endl;
  164. himag << ox << "\t" << oy << "\t" << depth << "\t"
  165. << std::imag(temp(0)) << "\t" << std::imag(temp(1))
  166. << "\t" << std::imag(temp(2)) << std::endl;
  167. depth += dx;
  168. }
  169. real.close();
  170. imag.close();
  171. hreal.close();
  172. himag.close();
  173. EmEarth->Delete();
  174. receivers->Delete();
  175. earth->Delete();
  176. //wire->Delete();
  177. #if LEMMAUSEVTK
  178. // Create the usual rendering stuff.
  179. vtkRenderer *renderer = vtkRenderer::New();
  180. vtkRenderWindow *renWin = vtkRenderWindow::New();
  181. vtkRenderWindowInteractor *iren = vtkRenderWindowInteractor::New();
  182. // Line of tx
  183. vtkLineSource *vline = vtkLineSource::New();
  184. vtkTubeFilter *vTube = vtkTubeFilter::New();
  185. vtkPolyDataMapper *vMapper = vtkPolyDataMapper::New();
  186. vtkActor *vActor = vtkActor::New();
  187. vline->SetPoint1(0,0,0);
  188. vline->SetPoint2(10,0,0);
  189. vTube->SetInputConnection(vline->GetOutputPort());
  190. vTube->SetRadius(.2);
  191. vTube->SetNumberOfSides(6);
  192. vMapper->SetInputConnection(vTube->GetOutputPort());
  193. vActor->SetMapper(vMapper);
  194. vActor->GetProperty()->SetColor(0.0, .0, 1.0);
  195. vActor->GetProperty()->SetOpacity(.15);
  196. renderer->AddActor(vActor);
  197. vtkLineSource *vline2 = vtkLineSource::New();
  198. vtkTubeFilter *vTube2 = vtkTubeFilter::New();
  199. vtkPolyDataMapper *vMapper2 = vtkPolyDataMapper::New();
  200. vtkActor *vActor2 = vtkActor::New();
  201. vline2->SetPoint1(10,0,0);
  202. vline2->SetPoint2(10,10,0);
  203. vTube2->SetInputConnection(vline2->GetOutputPort());
  204. vTube2->SetRadius(.2);
  205. vTube2->SetNumberOfSides(6);
  206. vMapper2->SetInputConnection(vTube2->GetOutputPort());
  207. vActor2->SetMapper(vMapper2);
  208. vActor2->GetProperty()->SetColor(0.0, .0, 1.0);
  209. vActor2->GetProperty()->SetOpacity(.15);
  210. renderer->AddActor(vActor2);
  211. vtkLineSource *vline3 = vtkLineSource::New();
  212. vtkTubeFilter *vTube3 = vtkTubeFilter::New();
  213. vtkPolyDataMapper *vMapper3 = vtkPolyDataMapper::New();
  214. vtkActor *vActor3 = vtkActor::New();
  215. vline3->SetPoint1(10,10,0);
  216. vline3->SetPoint2(0,10,0);
  217. vTube3->SetInputConnection(vline3->GetOutputPort());
  218. vTube3->SetRadius(.2);
  219. vTube3->SetNumberOfSides(6);
  220. vMapper3->SetInputConnection(vTube3->GetOutputPort());
  221. vActor3->SetMapper(vMapper3);
  222. vActor3->GetProperty()->SetColor(0.0, .0, 1.0);
  223. vActor3->GetProperty()->SetOpacity(.15);
  224. renderer->AddActor(vActor3);
  225. vtkLineSource *vline4 = vtkLineSource::New();
  226. vtkTubeFilter *vTube4 = vtkTubeFilter::New();
  227. vtkPolyDataMapper *vMapper4 = vtkPolyDataMapper::New();
  228. vtkActor *vActor4 = vtkActor::New();
  229. vline4->SetPoint1(0,10,0);
  230. vline4->SetPoint2(0,0,0);
  231. vTube4->SetInputConnection(vline4->GetOutputPort());
  232. vTube4->SetRadius(.2);
  233. vTube4->SetNumberOfSides(6);
  234. vMapper4->SetInputConnection(vTube4->GetOutputPort());
  235. vActor4->SetMapper(vMapper4);
  236. vActor4->GetProperty()->SetColor(0.0, .0, 1.0);
  237. vActor4->GetProperty()->SetOpacity(.15);
  238. renderer->AddActor(vActor4);
  239. loc << 50, 50, -1e-3;
  240. pa->ApproximateWithElectricDipoles(loc);
  241. vtkActor **pdipActors = new vtkActor*[pa->GetNumberOfDipoles()];
  242. std::cout << "Wire approximated with " << pa->GetNumberOfDipoles() << std::endl;
  243. for (int id=0; id<pa->GetNumberOfDipoles(); ++id) {
  244. pdipActors[id] = pa->GetVtkActor(id);
  245. renderer->AddActor(pdipActors[id]);
  246. }
  247. /*
  248. vtkActor **dipActors = new vtkActor*[wire->GetNumberOfDipoles()];
  249. for (int id=0; id<wire->GetNumberOfDipoles(); ++id) {
  250. dipActors[id] = wire->GetVtkActor(id);
  251. renderer->AddActor(dipActors[id]);
  252. }
  253. */
  254. renderer->SetBackground(1,1,1);
  255. // Render the window
  256. renWin->AddRenderer(renderer);
  257. renWin->SetWindowName("Wire antennae");
  258. iren->SetRenderWindow(renWin);
  259. iren->Initialize();
  260. iren->Start();
  261. iren->Render();
  262. #if 0
  263. cout << "Enter File name?: ";
  264. std::string pngName;
  265. std::cin >> pngName;
  266. vtkPNGWriter *pngwrite = vtkPNGWriter::New();
  267. vtkRenderLargeImage *renlarge = vtkRenderLargeImage::New();
  268. renlarge->SetInput(renderer);
  269. renlarge->SetMagnification(2);
  270. pngwrite->SetInputConnection(renlarge->GetOutputPort());
  271. pngName.append(".png");
  272. pngwrite->SetFileName(pngName.c_str());
  273. pngwrite->Write();
  274. #endif
  275. #endif // ----- not LEMMA_USE_VTK -----
  276. //std::cout << *pa << std::endl;
  277. //pa->Delete();
  278. return 0;
  279. }