Lemma is an Electromagnetics API
Você não pode selecionar mais de 25 tópicos Os tópicos devem começar com uma letra ou um número, podem incluir traços ('-') e podem ter até 35 caracteres.

utORS.cpp 10KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342
  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. #include "vtkRectilinearGrid.h"
  48. #include "vtkXMLRectilinearGridWriter.h"
  49. #include "vtkDoubleArray.h"
  50. #include "vtkFieldData.h"
  51. #include "vtkCellData.h"
  52. #endif // ----- not LEMMA_USE_VTK -----
  53. // For testing purposes disable VTK and run scale.sh
  54. //#undef LEMMAUSEVTK
  55. #include "timer.h"
  56. using namespace Lemma;
  57. double randDouble(double low, double high) {
  58. //srand(time(0));
  59. double temp;
  60. /* swap low & high around if the user makes no sense */
  61. if (low > high) {
  62. temp = low;
  63. low = high;
  64. high = temp;
  65. }
  66. /* calculate the random number & return it */
  67. temp = (rand() / (static_cast<double>(RAND_MAX) + 1.0))
  68. * (high - low) + low;
  69. return temp;
  70. }
  71. int main() {
  72. // Keep track of time
  73. jsw_timer timer;
  74. srand(time(0));
  75. PolygonalWireAntenna *pa = PolygonalWireAntenna::New();
  76. pa->SetNumberOfFrequencies(1);
  77. pa->SetFrequency(0, 540000);
  78. Real Depth = 370; // nominal depth
  79. Real Width = 0.14; // transmitter height
  80. Real Height = 0.70;// .14; //.014; // transmitter width
  81. pa->SetNumberOfPoints(5);
  82. /*
  83. pa->SetPoint(0, Vector3r( 0, 0, Depth));
  84. pa->SetPoint(1, Vector3r( 0, 0, Depth+Height));
  85. pa->SetPoint(2, Vector3r( 0, Width, Depth+Height));
  86. pa->SetPoint(3, Vector3r( 0, Width, Depth));
  87. pa->SetPoint(4, Vector3r( 0, 0, Depth));
  88. */
  89. pa->SetPoint(0, Vector3r( -Width/2., -Height/2., Depth));
  90. pa->SetPoint(1, Vector3r( Width/2., -Height/2., Depth));
  91. pa->SetPoint(2, Vector3r( Width/2., Height/2., Depth));
  92. pa->SetPoint(3, Vector3r( -Width/2., Height/2., Depth));
  93. pa->SetPoint(4, Vector3r( -Width/2., -Height/2., Depth));
  94. pa->SetCurrent(1.);
  95. pa->SetNumberOfTurns(6);
  96. //Vector3r rp = Vector3r::Random(3);
  97. //rp << 150., 10, 0.;
  98. //rp << -27.1456, 15.2350, -1e-3;
  99. //rp << randDouble(-35,35), randDouble(-35,35), randDouble(-35,35);
  100. //rp << 3.22806, -13.1548, 14.9695;
  101. //rp.setRandom(3);
  102. //std::cout << "rp " << rp.transpose() << std::endl;
  103. //pa->ApproximateWithElectricDipoles(rp);
  104. // Define model
  105. Real Sigma = 1.0/20. ; // .05;
  106. VectorXcr sigma(2);
  107. sigma << Complex(0.,0), Complex(Sigma,0);
  108. VectorXr thick(1);
  109. thick << 10;
  110. LayeredEarthEM *earth = LayeredEarthEM::New();
  111. earth->SetNumberOfLayers(2);
  112. earth->SetLayerConductivity(sigma);
  113. //earth->SetLayerThickness(thick);
  114. // Receivers
  115. ReceiverPoints *receivers = ReceiverPoints::New();
  116. Vector3r loc;
  117. Real ox = -5.*Width - .003373;
  118. Real oy = -3.*Height - .003373;
  119. Real oz = Depth - Height/2. - .003373;
  120. int nx = 120; // 60;
  121. int ny = 180; // 50;
  122. int nz = 100; // 40;
  123. Real hx = 11.*Width/nx;
  124. Real hy = 7.*Height/ny;
  125. Real hz = 1.*Height/nz;
  126. receivers->SetNumberOfReceivers(nx*ny*nz);
  127. int ir = 0;
  128. for (int iz=0; iz<nz; ++iz) {
  129. for (int iy=0; iy<ny; ++iy) {
  130. for (int ix=0; ix<nx; ++ix) {
  131. loc << ox+ix*hx, oy+iy*hy, oz+iz*hz;
  132. receivers->SetLocation(ir, loc);
  133. ++ ir;
  134. }
  135. }
  136. }
  137. // EmEarth
  138. EMEarth1D *EmEarth = EMEarth1D::New();
  139. //EmEarth->AttachWireAntenna(wire);
  140. EmEarth->AttachWireAntenna(pa);
  141. EmEarth->AttachLayeredEarthEM(earth);
  142. EmEarth->AttachReceiverPoints(receivers);
  143. EmEarth->SetFieldsToCalculate(E);
  144. //EmEarth->SetHankelTransformMethod(GAUSSIANQUADRATURE);
  145. EmEarth->SetHankelTransformMethod(ANDERSON801);
  146. // Do calculation
  147. timer.begin();
  148. EmEarth->CalculateWireAntennaFields();
  149. Real paTime = timer.end();
  150. std::cout << "Polygonal wire antennae time: " << paTime << "\n";
  151. //EmEarth->AttachWireAntenna(wire);
  152. //timer.begin();
  153. //EmEarth->CalculateWireAntennaFields();
  154. //Real waTime = timer.end();
  155. //std::cout << "Fixed wire antennae time: " << waTime << "\n";
  156. /*
  157. depth = depth2;
  158. std::fstream real("reale_lay.dat", std::ios::out);
  159. std::fstream imag("image_lay.dat", std::ios::out);
  160. for (int iz=0; iz<nz; ++iz) {
  161. Vector3cr temp = receivers->GetEfield(0,iz);
  162. real << ox << "\t" << oy << "\t" << depth << "\t"
  163. << temp(0).real() << "\t" << temp(1).real()
  164. << "\t" << temp(2).real() << std::endl;
  165. imag << ox << "\t" << oy << "\t" << depth << "\t"
  166. << std::imag(temp(0)) << "\t" << std::imag(temp(1))
  167. << "\t" << std::imag(temp(2)) << std::endl;
  168. depth += dx;
  169. }
  170. real.close();
  171. imag.close();
  172. */
  173. //wire->Delete();
  174. #if LEMMAUSEVTK
  175. // Set Coordinates
  176. vtkDoubleArray *xCoords = vtkDoubleArray::New();
  177. xCoords->InsertNextValue(ox-hx/2.);
  178. double xm1 = ox-hx/2.;
  179. for (int ix=0; ix<nx; ix++) {
  180. xCoords->InsertNextValue(xm1 + hx);
  181. xm1 += hx;
  182. }
  183. vtkDoubleArray *yCoords = vtkDoubleArray::New();
  184. yCoords->InsertNextValue(oy-hy/2.);
  185. double ym1 = oy-hy/2.;
  186. for (int iy=0; iy<ny; iy++) {
  187. yCoords->InsertNextValue(ym1 + hy);
  188. ym1 += hy;
  189. }
  190. vtkDoubleArray *zCoords = vtkDoubleArray::New();
  191. zCoords->InsertNextValue(oz-hz/2.);
  192. double zm1 = oz-hz/2.;
  193. for (int iz=0; iz<nz; iz++) {
  194. zCoords->InsertNextValue(zm1 + hz);
  195. zm1 += hz;
  196. }
  197. vtkDoubleArray *EReal = vtkDoubleArray::New();
  198. vtkDoubleArray *EImag = vtkDoubleArray::New();
  199. vtkDoubleArray *Watts = vtkDoubleArray::New();
  200. EReal->SetNumberOfComponents(3);
  201. EImag->SetNumberOfComponents(3);
  202. Watts->SetNumberOfComponents(1);
  203. ir = 0;
  204. Real WattsTotal(0);
  205. Real WattsInterior(0);
  206. for (int iz=0; iz<nz; ++iz) {
  207. for (int iy=0; iy<ny; ++iy) {
  208. for (int ix=0; ix<nx; ++ix) {
  209. //sigmaArray->InsertTuple1(i, sigma[ix][iy][iz] );
  210. Vector3cr E = receivers->GetEfield(0, ir);
  211. EReal-> InsertTuple3(ir, real(E(0)), real(E(1)), real(E(2)));
  212. EImag-> InsertTuple3(ir, imag(E(0)), imag(E(1)), imag(E(2)));
  213. //std::cout << std::abs(ox+ix*hx) << "\t" << Width/2. << endl; // && std::abs(oy+iy*hy) > Height/2. ) { // && std::abs(oz+iz*hz - Depth) > Width/2. ) {
  214. //if ( std::abs( ox+ix*hx ) < Width/2. && std::abs(oy+iy*hy) < Height/2. && std::abs(oz+iz*hz - Depth) < Width/2. ) {
  215. if ( std::sqrt( std::pow(ox+ix*hx,2) + std::pow(oz+iz*hz - Depth, 2) ) < .085 && std::abs(oy+iy*hy) < Height/2. ) {
  216. Watts-> InsertTuple1(ir, 1e-20 );
  217. WattsInterior += .5* (( pow((std::abs(E(0)) + std::abs(E(1)) + std::abs(E(2))), 2)*Sigma)*hx*hy*hz );
  218. } else {
  219. Watts-> InsertTuple1(ir, .5* (pow( (std::abs(E(0)) + std::abs(E(1)) + std::abs(E(2))), 2)*Sigma)*hx*hy*hz );
  220. //Watts-> InsertTuple1(ir, 1e-20 );
  221. WattsTotal += .5* (( pow((std::abs(E(0)) + std::abs(E(1)) + std::abs(E(2))), 2)*Sigma)*hx*hy*hz );
  222. }
  223. ++ ir;
  224. }
  225. }
  226. }
  227. std::cout << "Total Power: " << Sigma << "\t" << WattsTotal << "\t" << WattsInterior << endl;
  228. EReal->SetName("E_real");
  229. EImag->SetName("E_imag");
  230. Watts->SetName("Power");
  231. vtkRectilinearGrid *rgrid = vtkRectilinearGrid::New();
  232. rgrid->SetDimensions(nx+1,ny+1,nz+1);
  233. rgrid->SetXCoordinates(xCoords);
  234. rgrid->SetYCoordinates(yCoords);
  235. rgrid->SetZCoordinates(zCoords);
  236. rgrid->GetCellData()->AddArray(EReal);
  237. rgrid->GetCellData()->AddArray(EImag);
  238. rgrid->GetCellData()->AddArray(Watts);
  239. //rgrid->Update();
  240. vtkXMLRectilinearGridWriter *gridWrite = vtkXMLRectilinearGridWriter::New();
  241. gridWrite->SetInputData(rgrid);
  242. gridWrite->SetFileName("ors.vtr");
  243. gridWrite->Write();
  244. //gridWrite->Update();
  245. #endif
  246. EmEarth->Delete();
  247. receivers->Delete();
  248. earth->Delete();
  249. #if LEMMAUSEVTKX
  250. // Create the usual rendering stuff.
  251. vtkRenderer *renderer = vtkRenderer::New();
  252. vtkRenderWindow *renWin = vtkRenderWindow::New();
  253. vtkRenderWindowInteractor *iren = vtkRenderWindowInteractor::New();
  254. std::cout << "wire antennae approximating " << std::endl;
  255. loc << 0, .5*Width, Depth+.5*Height;
  256. pa->ApproximateWithElectricDipoles(loc);
  257. std::cout << "Wire approximated with " << pa->GetNumberOfDipoles() << std::endl;
  258. vtkActor **pdipActors = new vtkActor*[pa->GetNumberOfDipoles()];
  259. for (int id=0; id<pa->GetNumberOfDipoles(); ++id) {
  260. pdipActors[id] = pa->GetVtkActor(id);
  261. renderer->AddActor(pdipActors[id]);
  262. }
  263. /*
  264. vtkActor **dipActors = new vtkActor*[wire->GetNumberOfDipoles()];
  265. for (int id=0; id<wire->GetNumberOfDipoles(); ++id) {
  266. dipActors[id] = wire->GetVtkActor(id);
  267. renderer->AddActor(dipActors[id]);
  268. }
  269. */
  270. renderer->SetBackground(1,1,1);
  271. // Render the window
  272. renWin->AddRenderer(renderer);
  273. renWin->SetWindowName("Wire antennae");
  274. iren->SetRenderWindow(renWin);
  275. iren->Initialize();
  276. iren->Start();
  277. iren->Render();
  278. #if 0
  279. cout << "Enter File name?: ";
  280. std::string pngName;
  281. std::cin >> pngName;
  282. vtkPNGWriter *pngwrite = vtkPNGWriter::New();
  283. vtkRenderLargeImage *renlarge = vtkRenderLargeImage::New();
  284. renlarge->SetInput(renderer);
  285. renlarge->SetMagnification(2);
  286. pngwrite->SetInputConnection(renlarge->GetOutputPort());
  287. pngName.append(".png");
  288. pngwrite->SetFileName(pngName.c_str());
  289. pngwrite->Write();
  290. #endif
  291. #endif // ----- not LEMMA_USE_VTK -----
  292. //std::cout << *pa << std::endl;
  293. //pa->Delete();
  294. return 0;
  295. }