3D EM based on Schur decomposition
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.

EMSchur3DBase.cpp 59KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518
  1. // ===========================================================================
  2. //
  3. // Filename: EMSchur3DBase.cpp
  4. //
  5. // Created: 09/20/2013 04:53:40 PM
  6. // Compiler: Tested with g++, icpc, and MSVC 2010
  7. //
  8. // Author: Trevor Irons (ti)
  9. //
  10. // Organisation: University of Utah
  11. //
  12. // Email: tirons@egi.utah.edu
  13. //
  14. // ===========================================================================
  15. /**
  16. @file
  17. @author Trevor Irons
  18. @date 09/20/2013
  19. @version $Id$
  20. THis is a port of the procedural code developed at School of Mines
  21. **/
  22. #include "EMSchur3DBase.h"
  23. typedef Eigen::Triplet<Lemma::Complex> Tc;
  24. typedef Eigen::Triplet<Lemma::Real> Tr;
  25. #define UPPER 0 // LOWER WAS 0
  26. #define LOWER 1 // 1=true, 0=false
  27. namespace Lemma {
  28. // ==================== FRIEND METHODS =====================
  29. std::ostream &operator<<(std::ostream &stream, const EMSchur3DBase &ob) {
  30. stream << ob.Serialize() << "\n";
  31. return stream;
  32. }
  33. // ==================== LIFECYCLE =======================
  34. //--------------------------------------------------------------------------------------
  35. // Class: EMSchur3DBase
  36. // Method: EMSchur3DBase
  37. // Description: constructor (protected)
  38. //--------------------------------------------------------------------------------------
  39. // std::shared_ptr<EMSchur3DBase> EMSchur3DBase::NewSP() {
  40. // return std::make_shared<EMSchur3DBase>( ctor_key() );
  41. // }
  42. //--------------------------------------------------------------------------------------
  43. // Class: EMSchur3DBase
  44. // Method: EMSchur3DBase
  45. // Description: constructor (protected)
  46. //--------------------------------------------------------------------------------------
  47. EMSchur3DBase::EMSchur3DBase ( const ctor_key& key ) : LemmaObject( key ),
  48. Grid(nullptr),
  49. //Survey(nullptr),
  50. LayModel(nullptr), Cvec(nullptr),
  51. ResFile("source"), sigma(nullptr), sigmap(nullptr)
  52. {
  53. } // ----- end of method EMSchur3DBase::EMSchur3DBase (constructor) -----
  54. //--------------------------------------------------------------------------------------
  55. // Class: EMSchur3DBase
  56. // Method: EMSchur3DBase
  57. // Description: constructor (protected)
  58. //--------------------------------------------------------------------------------------
  59. EMSchur3DBase::EMSchur3DBase ( const YAML::Node& node, const ctor_key& key ) : LemmaObject( key ),
  60. Grid(nullptr),
  61. //Survey(nullptr),
  62. LayModel(nullptr), Cvec(nullptr),
  63. ResFile("source"), sigma(nullptr), sigmap(nullptr)
  64. {
  65. } // ----- end of method EMSchur3DBase::~EMSchur3DBase (destructor) -----
  66. //--------------------------------------------------------------------------------------
  67. // Class: EMSchur3DBase
  68. // Method: ~EMSchur3DBase
  69. // Description: destructor (protected)
  70. //--------------------------------------------------------------------------------------
  71. EMSchur3DBase::~EMSchur3DBase ( ) {
  72. if (sigma) Delete3DScalar(sigma);
  73. if (sigmap) Delete3DScalar(sigmap);
  74. if (Cvec) delete [] Cvec;
  75. //if (CvecRe) delete [] CvecRe;
  76. } // ----- end of method EMSchur3DBase::~EMSchur3DBase (destructor) -----
  77. //--------------------------------------------------------------------------------------
  78. // Class: EMSchur3D
  79. // Method: Serialize
  80. //--------------------------------------------------------------------------------------
  81. YAML::Node EMSchur3DBase::Serialize ( ) const {
  82. YAML::Node node = LemmaObject::Serialize();
  83. node.SetTag( this->GetName() );
  84. node["EMSchur3D_VERSION"] = EMSCHUR3D_VERSION;
  85. return node;
  86. } // ----- end of method EMSchur3D::Serialize -----
  87. //--------------------------------------------------------------------------------------
  88. // Class: EMSchur3DBase
  89. // Method: BuildC
  90. //--------------------------------------------------------------------------------------
  91. void EMSchur3DBase::BuildC ( Real*** sigmax, Real*** sigmay, Real*** sigmaz, const int& iw) {
  92. Cvec[iw].resize( unx+uny+unz , unx+uny+unz );
  93. #if LOWER && UPPER
  94. Cvec[iw].reserve(Eigen::VectorXi::Constant(unx+uny+unz, 7)); // Whole
  95. #else
  96. Cvec[iw].reserve(Eigen::VectorXi::Constant(unx+uny+unz, 4)); // Upper/Lower
  97. #endif
  98. //Cvec_s.resize( idx. )
  99. //CMMvec[iw].resize( unx+uny+unz , unx+uny+unz );
  100. //CMMvec[iw].reserve(Eigen::VectorXi::Constant(unx+uny+unz, 1));
  101. Real omega = Omegas[iw];
  102. std::cout << "Building C_" << iw << std::endl;
  103. //Complex hsig(0);
  104. Real hsig(0);
  105. Real EPS(1e-24);
  106. Real EPS2(1e-18);
  107. // book keeping
  108. int ir = 0;
  109. int ic = 0;
  110. // LAPL{Ax} + i omega mu sigmax
  111. for (int iz=0; iz<nz; ++iz) {
  112. for (int iy=0; iy<ny; ++iy) {
  113. for (int ix=0; ix<nx+1; ++ix) {
  114. // Calculate 1/2 step values on staggered grid
  115. Real alo_ihz2(0);
  116. Real ahi_ihz2(0);
  117. if (iz == 0) {
  118. ahi_ihz2 = std::pow(.5*hz[iz+1] + .5*hz[iz], -2);
  119. alo_ihz2 = std::pow(.5*hz[iz+1] + .5*hz[iz], -2);
  120. } else if (iz == nz-1) {
  121. alo_ihz2 = std::pow(.5*hz[iz-1] + .5*hz[iz], -2);
  122. ahi_ihz2 = std::pow(.5*hz[iz-1] + .5*hz[iz], -2);
  123. } else {
  124. alo_ihz2 = std::pow(.5*hz[iz-1] + .5*hz[iz], -2);
  125. ahi_ihz2 = std::pow(.5*hz[iz+1] + .5*hz[iz], -2);
  126. }
  127. Real alo_ihy2(0); // half step low jump in y
  128. Real ahi_ihy2(0); // half step high jump in y
  129. if (iy == 0) {
  130. ahi_ihy2 = std::pow(.5*hy[iy+1] + .5*hy[iy], -2);
  131. alo_ihy2 = std::pow(.5*hy[iy+1] + .5*hy[iy], -2);
  132. } else if (iy == ny-1) {
  133. alo_ihy2 = std::pow(.5*hy[iy-1] + .5*hy[iy], -2);
  134. ahi_ihy2 = std::pow(.5*hy[iy-1] + .5*hy[iy], -2);
  135. } else {
  136. alo_ihy2 = std::pow(.5*hy[iy-1] + .5*hy[iy], -2);
  137. ahi_ihy2 = std::pow(.5*hy[iy+1] + .5*hy[iy], -2);
  138. }
  139. /////////////////////////////////////////////////////////////////////////
  140. // Diagonal entry
  141. // Harmonically average sigmax
  142. if (ix == 0) {
  143. hsig = sigma[ix][iy][iz];
  144. } else if (ix == nx) {
  145. hsig = sigma[ix-1][iy][iz];
  146. } else {
  147. hsig = ((hx[ix]+hx[ix-1])/2.) *
  148. (1. / ( (hx[ix ] / (2.*sigma[ix ][iy][iz] + EPS)) +
  149. ( hx[ix-1] / (2.*sigma[ix-1][iy][iz] + EPS)) ) );
  150. }
  151. if (std::abs(hsig) < EPS2) hsig = 0;
  152. //hsig += Complex(0, omega*EPSILON0);
  153. Real Sum(0);
  154. if (ix == 0) {
  155. // Dirichlet bdry condition on A
  156. // Sum = alo_ihy2+ahi_ihy2 + alo_ihz2+ahi_ihz2 + 2.*ihx2[ix ];
  157. // Neumann bdry
  158. Sum = alo_ihy2+ahi_ihy2 + alo_ihz2+ahi_ihz2 + ihx2[ix ];
  159. } else if (ix == nx) {
  160. // Dirichlet bdry condition on A
  161. // Sum = alo_ihy2+ahi_ihy2 + alo_ihz2+ahi_ihz2 + 2.*ihx2[ix-1];
  162. // Neumann bdry
  163. Sum = alo_ihy2+ahi_ihy2 + alo_ihz2+ahi_ihz2 + ihx2[ix-1];
  164. } else {
  165. Sum = alo_ihy2+ahi_ihy2 + alo_ihz2+ahi_ihz2 + ihx2[ix] + ihx2[ix-1];
  166. }
  167. /////////////////////////////////////////
  168. // minus side off diagonal terms
  169. #if LOWER
  170. // Third Off Diagonal
  171. if (iz!=0) Cvec[iw].insert(ir, ic-ny*(nx+1)) = Complex(-ahi_ihz2, 0);
  172. // Second Off Diagonal
  173. if (iy!=0) Cvec[iw].insert(ir, ic-(nx+1)) = Complex(-ahi_ihy2, 0);
  174. // First Off Diagonal
  175. if (ix!=0) Cvec[iw].insert(ir, ic-1) = Complex(-ihx2[ix-1], 0);
  176. #endif
  177. ////////////////////////////////////////////
  178. // Diagonal Term
  179. Cvec[iw].insert(ir,ic) = Complex(Sum, 0) + Complex(0,1)*omega*MU0*hsig;
  180. //Cvec[iw].insert(ir,ic) = Complex(Sum, omega*MU0*hsig);
  181. //CMMvec[iw].insert(ir,ic) = 1./ Complex(Sum, omega*MU0*hsig);
  182. ////////////////////////////////////////////////////////////////////////
  183. // plus side off diagonal terms
  184. #if UPPER
  185. // First Off Diagonal
  186. if (ix!=nx) Cvec[iw].insert(ir, ic+1) = Complex(-ihx2[ix], 0);
  187. // Second Off Diagonal
  188. if (iy!=ny-1) Cvec[iw].insert(ir, ic+(nx+1)) = Complex(-ahi_ihy2, 0);
  189. // Third Off Diagonal
  190. if (iz!=nz-1) Cvec[iw].insert(ir, ic+ny*(nx+1)) = Complex(-ahi_ihz2, 0);
  191. #endif
  192. ++ir;
  193. ++ic;
  194. }
  195. }
  196. }
  197. assert(ic == unx);
  198. // LAPL{Ay} + i omega mu sigmay
  199. for (int iz=0; iz<nz; ++iz) {
  200. for (int iy=0; iy<ny+1; ++iy) {
  201. for (int ix=0; ix<nx; ++ix) {
  202. // Calculate 1/2 step values on staggered grid
  203. Real alo_ihz2(0);
  204. Real ahi_ihz2(0);
  205. if (iz == 0) {
  206. ahi_ihz2 = std::pow(.5*hz[iz+1] + .5*hz[iz], -2);
  207. alo_ihz2 = std::pow(.5*hz[iz+1] + .5*hz[iz], -2);
  208. } else if (iz == nz-1) {
  209. alo_ihz2 = std::pow(.5*hz[iz-1] + .5*hz[iz], -2);
  210. ahi_ihz2 = std::pow(.5*hz[iz-1] + .5*hz[iz], -2);
  211. } else {
  212. alo_ihz2 = std::pow(.5*hz[iz-1] + .5*hz[iz], -2);
  213. ahi_ihz2 = std::pow(.5*hz[iz+1] + .5*hz[iz], -2);
  214. }
  215. // Calculate 1/2 step values on staggered grid
  216. Real alo_ihx2(0);
  217. Real ahi_ihx2(0);
  218. if (ix == 0) {
  219. ahi_ihx2 = std::pow(.5*hx[ix+1] + .5*hx[ix], -2);
  220. alo_ihx2 = std::pow(.5*hx[ix+1] + .5*hx[ix], -2);
  221. } else if (ix == nx-1) {
  222. alo_ihx2 = std::pow(.5*hx[ix-1] + .5*hx[ix], -2);
  223. ahi_ihx2 = std::pow(.5*hx[ix-1] + .5*hx[ix], -2);
  224. } else {
  225. alo_ihx2 = std::pow(.5*hx[ix-1] + .5*hx[ix], -2);
  226. ahi_ihx2 = std::pow(.5*hx[ix+1] + .5*hx[ix], -2);
  227. }
  228. // Harmonically average sigmay
  229. if (iy == 0) {
  230. hsig = sigma[ix][iy][iz];
  231. } else if (iy == ny) {
  232. hsig = sigma[ix][iy-1][iz];
  233. } else {
  234. hsig = ((hy[iy]+hy[iy-1])/2.) *
  235. (1. / ( (hy[iy ] / (2.*sigma[ix][iy ][iz] + EPS)) +
  236. ( hy[iy-1] / (2.*sigma[ix][iy-1][iz] + EPS)) ) );
  237. }
  238. if (std::abs(hsig) < EPS2) hsig = 0;
  239. //hsig += Complex(0, omega*EPSILON0);
  240. Real Sum(0);
  241. if (iy == 0) {
  242. // Dirichlet bdry condition on A
  243. //Sum = alo_ihx2+ahi_ihx2 + alo_ihz2+ahi_ihz2 + 2.*ihy2[iy ];
  244. // Neumann bdry
  245. Sum = alo_ihx2+ahi_ihx2 + alo_ihz2+ahi_ihz2 + ihy2[iy ];
  246. } else if (iy == ny) {
  247. // Dirichlet bdry condition on A
  248. //Sum = alo_ihx2+ahi_ihx2 + alo_ihz2+ahi_ihz2 + 2.*ihy2[iy-1];
  249. // Neumann bdry
  250. Sum = alo_ihx2+ahi_ihx2 + alo_ihz2+ahi_ihz2 + ihy2[iy-1];
  251. } else {
  252. Sum = alo_ihx2+ahi_ihx2 + alo_ihz2+ahi_ihz2 + ihy2[iy] + ihy2[iy-1];
  253. }
  254. #if LOWER
  255. // Third Off Diagonal
  256. if (iz!=0) Cvec[iw].insert(ir,ic-(ny+1)*(nx)) = Complex(-ahi_ihz2, 0);
  257. // Second Off Diagonal
  258. if (iy!=0) Cvec[iw].insert(ir,ic-(nx)) = Complex(-ihy2[iy-1], 0);
  259. // First Off Diagonal
  260. if (ix!=0) Cvec[iw].insert(ir,ic-1) = Complex(-ahi_ihx2, 0);
  261. #endif
  262. // Diagonal Term
  263. Cvec[iw].insert(ir,ic) = Complex(Sum, 0) + Complex(0,1)*omega*MU0*hsig;
  264. //Cvec[iw].insert(ir,ic) = Complex(Sum, omega*MU0*hsig);
  265. //CMMvec[iw].insert(ir,ic) = 1./ Complex(Sum, omega*MU0*hsig);
  266. #if UPPER
  267. // First Off Diagonal
  268. if (ix!=nx-1) Cvec[iw].insert(ir,ic+1) = Complex(-ahi_ihx2, 0);
  269. // Second Off Diagonal
  270. if (iy!=ny) Cvec[iw].insert(ir,ic+(nx)) = Complex(-ihy2[iy], 0);
  271. // Third Off Diagonal
  272. if (iz!=nz-1) Cvec[iw].insert(ir,ic+(ny+1)*(nx)) = Complex(-ahi_ihz2, 0);
  273. #endif
  274. ++ir;
  275. ++ic;
  276. }
  277. }
  278. }
  279. assert(ic == unx+uny);
  280. // LAPL{Az} + i omega mu sigmaz
  281. for (int iz=0; iz<nz+1; ++iz) {
  282. for (int iy=0; iy<ny; ++iy) {
  283. for (int ix=0; ix<nx; ++ix) {
  284. // Calculate 1/2 step values on staggered grid
  285. Real alo_ihx2(0);
  286. Real ahi_ihx2(0);
  287. if (ix == 0) {
  288. ahi_ihx2 = std::pow(.5*hx[ix+1] + .5*hx[ix], -2);
  289. alo_ihx2 = std::pow(.5*hx[ix+1] + .5*hx[ix], -2);
  290. } else if (ix == nx-1) {
  291. alo_ihx2 = std::pow(.5*hx[ix-1] + .5*hx[ix], -2);
  292. ahi_ihx2 = std::pow(.5*hx[ix-1] + .5*hx[ix], -2);
  293. } else {
  294. alo_ihx2 = std::pow(.5*hx[ix-1] + .5*hx[ix], -2);
  295. ahi_ihx2 = std::pow(.5*hx[ix+1] + .5*hx[ix], -2);
  296. }
  297. // Calculate 1/2 step values on staggered grid
  298. Real alo_ihy2(0);
  299. Real ahi_ihy2(0);
  300. if (iy == 0) {
  301. ahi_ihy2 = std::pow(.5*hy[iy+1] + .5*hy[iy], -2);
  302. alo_ihy2 = std::pow(.5*hy[iy+1] + .5*hy[iy], -2);
  303. } else if (iy == ny-1) {
  304. alo_ihy2 = std::pow(.5*hy[iy-1] + .5*hy[iy], -2);
  305. ahi_ihy2 = std::pow(.5*hy[iy-1] + .5*hy[iy], -2);
  306. } else {
  307. alo_ihy2 = std::pow(.5*hy[iy-1] + .5*hy[iy], -2);
  308. ahi_ihy2 = std::pow(.5*hy[iy+1] + .5*hy[iy], -2);
  309. }
  310. // Harmonically average sigmaz
  311. if (iz==0) {
  312. hsig = sigma[ix][iy][iz];
  313. } else if (iz == nz) {
  314. hsig = sigma[ix][iy][iz-1];
  315. } else {
  316. hsig = ((hz[iz]+hz[iz-1])/2.) *
  317. (1. / ( (hz[iz ] / (2.*sigma[ix][iy][iz ] + EPS)) +
  318. ( hz[iz-1] / (2.*sigma[ix][iy][iz-1] + EPS)) ) );
  319. }
  320. if (std::abs(hsig) < EPS2) hsig = 0;
  321. //hsig += Complex(0, omega*EPSILON0);
  322. Real Sum(0);
  323. if (iz == 0) {
  324. // Dirichlet bdry condition on A
  325. //Sum = alo_ihx2 + ahi_ihx2 + alo_ihy2 + ahi_ihy2 + 2.*ihz2[iz ];
  326. // Neumann bdry
  327. Sum = alo_ihx2 + ahi_ihx2 + alo_ihy2 + ahi_ihy2 + ihz2[iz ];
  328. } else if (iz == nz) {
  329. // Dirichlet bdry condition on A
  330. Sum = alo_ihx2 + ahi_ihx2 + alo_ihy2 + ahi_ihy2 + 2.*ihz2[iz-1];
  331. // Neumann bdry
  332. //Sum = alo_ihx2 + ahi_ihx2 + alo_ihy2 + ahi_ihy2 + ihz2[iz-1];
  333. } else {
  334. Sum = alo_ihx2 + ahi_ihx2 + alo_ihy2 + ahi_ihy2 + ihz2[iz] + ihz2[iz-1];
  335. }
  336. #if LOWER
  337. // Third Off Diagonal
  338. if (iz!=0) Cvec[iw].insert(ir,ic-ny*(nx)) = Complex(-ihz2[iz-1], 0);
  339. // Second Off Diagonal
  340. if (iy!=0) Cvec[iw].insert(ir,ic-(nx)) = Complex(-ahi_ihy2, 0);
  341. // First Off Diagonal
  342. if (ix!=0) Cvec[iw].insert(ir,ic-1) = Complex(-ahi_ihx2, 0);
  343. #endif
  344. // Diagonal Term
  345. Cvec[iw].insert(ir,ic) = Complex(Sum, 0) + Complex(0,1)*omega*MU0*hsig;
  346. //Cvec[iw].insert(ir,ic) = Complex(Sum, omega*MU0*hsig);
  347. //CMMvec[iw].insert(ir,ic) = 1. / Complex(Sum, omega*MU0*hsig);
  348. #if UPPER
  349. // First Off Diagonal
  350. if (ix!=nx-1) Cvec[iw].insert(ir,ic+1) = Complex(-ahi_ihx2, 0);
  351. // Second Off Diagonal
  352. if (iy!=ny-1) Cvec[iw].insert(ir,ic+(nx)) = Complex(-ahi_ihy2, 0);
  353. // Third Off Diagonal
  354. if (iz!=nz) Cvec[iw].insert(ir,ic+ny*(nx)) = Complex(-ihz2[iz], 0);
  355. #endif
  356. ++ir;
  357. ++ic;
  358. }
  359. }
  360. }
  361. assert(ic == unx+uny+unz);
  362. //Cvec[iw].makeCompressed();
  363. Cvec[iw].makeCompressed();
  364. //CMMvec[iw].makeCompressed();
  365. //CsymVec[iw] = Cvec[iw].selfadjointView<Eigen::Upper>();
  366. return ;
  367. } // ----- end of method EMSchur3DBase::BuildC -----
  368. //--------------------------------------------------------------------------------------
  369. // Class: EMSchur3DBase
  370. // Method: BuildCReal
  371. //--------------------------------------------------------------------------------------
  372. /*
  373. void EMSchur3DBase::BuildCReal ( Real*** sigmax, Real*** sigmay, Real*** sigmaz, const int& iw) {
  374. int CI = unx+uny+unz;
  375. CvecRe[iw].resize( 2*CI, 2*CI ); // twice the size, but Real
  376. std::vector<Tr> triplets;
  377. triplets.reserve( 9*CI );
  378. Real omega = Omegas[iw];
  379. std::cout << "Building real C_" << iw << std::endl;
  380. // LAPL{Ax} + i omega mu sigmax
  381. int ir = 0;
  382. int ic = 0;
  383. Real hsig(0);
  384. Real EPS(1e-24);
  385. Real EPS2(1e-18);
  386. for (int iz=0; iz<nz; ++iz) {
  387. for (int iy=0; iy<ny; ++iy) {
  388. for (int ix=0; ix<nx+1; ++ix) {
  389. // Calculate 1/2 step values on staggered grid
  390. Real alo_ihz2(0);
  391. Real ahi_ihz2(0);
  392. if (iz == 0) {
  393. ahi_ihz2 = std::pow(.5*hz[iz+1] + .5*hz[iz], -2);
  394. alo_ihz2 = std::pow(.5*hz[iz+1] + .5*hz[iz], -2);
  395. } else if (iz == nz-1) {
  396. alo_ihz2 = std::pow(.5*hz[iz-1] + .5*hz[iz], -2);
  397. ahi_ihz2 = std::pow(.5*hz[iz-1] + .5*hz[iz], -2);
  398. } else {
  399. alo_ihz2 = std::pow(.5*hz[iz-1] + .5*hz[iz], -2);
  400. ahi_ihz2 = std::pow(.5*hz[iz+1] + .5*hz[iz], -2);
  401. }
  402. Real alo_ihy2(0); // half step low jump in y
  403. Real ahi_ihy2(0); // half step high jump in y
  404. if (iy == 0) {
  405. ahi_ihy2 = std::pow(.5*hy[iy+1] + .5*hy[iy], -2);
  406. alo_ihy2 = std::pow(.5*hy[iy+1] + .5*hy[iy], -2);
  407. } else if (iy == ny-1) {
  408. alo_ihy2 = std::pow(.5*hy[iy-1] + .5*hy[iy], -2);
  409. ahi_ihy2 = std::pow(.5*hy[iy-1] + .5*hy[iy], -2);
  410. } else {
  411. alo_ihy2 = std::pow(.5*hy[iy-1] + .5*hy[iy], -2);
  412. ahi_ihy2 = std::pow(.5*hy[iy+1] + .5*hy[iy], -2);
  413. }
  414. /////////////////////////////////////////////////////////////////////////
  415. // Diagonal entry
  416. // Harmonically average sigmax
  417. if (ix == 0) {
  418. hsig = sigma[ix][iy][iz];
  419. } else if (ix == nx) {
  420. hsig = sigma[ix-1][iy][iz];
  421. } else {
  422. hsig = ((hx[ix]+hx[ix-1])/2.) *
  423. (1. / ( (hx[ix ] / (2.*sigma[ix ][iy][iz] + EPS)) +
  424. ( hx[ix-1] / (2.*sigma[ix-1][iy][iz] + EPS)) ) );
  425. }
  426. if (std::abs(hsig) < EPS2) hsig = 0;
  427. //hsig += Complex(0, omega*EPSILON0);
  428. Real Sum(0);
  429. if (ix == 0) {
  430. // Dirichlet bdry condition on A
  431. // Sum = alo_ihy2+ahi_ihy2 + alo_ihz2+ahi_ihz2 + 2.*ihx2[ix ];
  432. // Neumann bdry
  433. Sum = alo_ihy2+ahi_ihy2 + alo_ihz2+ahi_ihz2 + ihx2[ix ];
  434. } else if (ix == nx) {
  435. // Dirichlet bdry condition on A
  436. // Sum = alo_ihy2+ahi_ihy2 + alo_ihz2+ahi_ihz2 + 2.*ihx2[ix-1];
  437. // Neumann bdry
  438. Sum = alo_ihy2+ahi_ihy2 + alo_ihz2+ahi_ihz2 + ihx2[ix-1];
  439. } else {
  440. Sum = alo_ihy2+ahi_ihy2 + alo_ihz2+ahi_ihz2 + ihx2[ix] + ihx2[ix-1];
  441. }
  442. ////////////////////////////////////////////
  443. // Diagonal Term
  444. triplets.push_back( Tr(ir, ic , Sum) );
  445. triplets.push_back( Tr(ir+CI, ic+CI, -Sum) );
  446. ////////////////////////////////////////////////////////////////////////
  447. // plus side off diagonal terms
  448. // First Off Diagonal
  449. if (ix!=nx) {
  450. triplets.push_back( Tr(ir, ic+1 , -ihx2[ix]) );
  451. triplets.push_back( Tr(ir+CI, ic+1+CI, ihx2[ix]) );
  452. }
  453. // Second Off Diagonal
  454. if (iy!=ny-1) {
  455. triplets.push_back( Tr(ir , ic+(nx+1) , -ahi_ihy2) );
  456. triplets.push_back( Tr(ir+CI, ic+(nx+1)+CI, ahi_ihy2) );
  457. }
  458. // Third Off Diagonal
  459. if (iz!=nz-1) {
  460. triplets.push_back( Tr(ir , ic+ny*(nx+1) , -ahi_ihz2) );
  461. triplets.push_back( Tr(ir+CI, ic+ny*(nx+1)+CI, ahi_ihz2) );
  462. }
  463. ////////////////////////////////////////////////////////////////////////
  464. // imaginary part
  465. triplets.push_back( Tr(ir, ic+CI, omega*MU0*hsig) );
  466. ++ir;
  467. ++ic;
  468. }
  469. }
  470. }
  471. assert(ic == unx);
  472. // LAPL{Ay} + i omega mu sigmay
  473. for (int iz=0; iz<nz; ++iz) {
  474. for (int iy=0; iy<ny+1; ++iy) {
  475. for (int ix=0; ix<nx; ++ix) {
  476. // Calculate 1/2 step values on staggered grid
  477. Real alo_ihz2(0);
  478. Real ahi_ihz2(0);
  479. if (iz == 0) {
  480. ahi_ihz2 = std::pow(.5*hz[iz+1] + .5*hz[iz], -2);
  481. alo_ihz2 = std::pow(.5*hz[iz+1] + .5*hz[iz], -2);
  482. } else if (iz == nz-1) {
  483. alo_ihz2 = std::pow(.5*hz[iz-1] + .5*hz[iz], -2);
  484. ahi_ihz2 = std::pow(.5*hz[iz-1] + .5*hz[iz], -2);
  485. } else {
  486. alo_ihz2 = std::pow(.5*hz[iz-1] + .5*hz[iz], -2);
  487. ahi_ihz2 = std::pow(.5*hz[iz+1] + .5*hz[iz], -2);
  488. }
  489. // Calculate 1/2 step values on staggered grid
  490. Real alo_ihx2(0);
  491. Real ahi_ihx2(0);
  492. if (ix == 0) {
  493. ahi_ihx2 = std::pow(.5*hx[ix+1] + .5*hx[ix], -2);
  494. alo_ihx2 = std::pow(.5*hx[ix+1] + .5*hx[ix], -2);
  495. } else if (ix == nx-1) {
  496. alo_ihx2 = std::pow(.5*hx[ix-1] + .5*hx[ix], -2);
  497. ahi_ihx2 = std::pow(.5*hx[ix-1] + .5*hx[ix], -2);
  498. } else {
  499. alo_ihx2 = std::pow(.5*hx[ix-1] + .5*hx[ix], -2);
  500. ahi_ihx2 = std::pow(.5*hx[ix+1] + .5*hx[ix], -2);
  501. }
  502. // Harmonically average sigmay
  503. if (iy == 0) {
  504. hsig = sigma[ix][iy][iz];
  505. } else if (iy == ny) {
  506. hsig = sigma[ix][iy-1][iz];
  507. } else {
  508. hsig = ((hy[iy]+hy[iy-1])/2.) *
  509. (1. / ( (hy[iy ] / (2.*sigma[ix][iy ][iz] + EPS)) +
  510. ( hy[iy-1] / (2.*sigma[ix][iy-1][iz] + EPS)) ) );
  511. }
  512. if (std::abs(hsig) < EPS2) hsig = 0;
  513. //hsig += Complex(0, omega*EPSILON0);
  514. Real Sum(0);
  515. if (iy == 0) {
  516. // Dirichlet bdry condition on A
  517. //Sum = alo_ihx2+ahi_ihx2 + alo_ihz2+ahi_ihz2 + 2.*ihy2[iy ];
  518. // Neumann bdry
  519. Sum = alo_ihx2+ahi_ihx2 + alo_ihz2+ahi_ihz2 + ihy2[iy ];
  520. } else if (iy == ny) {
  521. // Dirichlet bdry condition on A
  522. //Sum = alo_ihx2+ahi_ihx2 + alo_ihz2+ahi_ihz2 + 2.*ihy2[iy-1];
  523. // Neumann bdry
  524. Sum = alo_ihx2+ahi_ihx2 + alo_ihz2+ahi_ihz2 + ihy2[iy-1];
  525. } else {
  526. Sum = alo_ihx2+ahi_ihx2 + alo_ihz2+ahi_ihz2 + ihy2[iy] + ihy2[iy-1];
  527. }
  528. //////////////////////////////////////////////////////////////////////////////
  529. // Diagonal Term
  530. triplets.push_back( Tr(ir, ic , Sum) );
  531. triplets.push_back( Tr(ir+CI, ic+CI, -Sum) );
  532. // First Off Diagonal
  533. if (ix!=nx-1) {
  534. triplets.push_back( Tr(ir, ic+1 , -ahi_ihx2));
  535. triplets.push_back( Tr(ir+CI, ic+1+CI, ahi_ihx2));
  536. }
  537. // Second Off Diagonal
  538. if (iy!=ny) {
  539. triplets.push_back( Tr(ir , ic+(nx) , -ihy2[iy]) );
  540. triplets.push_back( Tr(ir+CI, ic+(nx)+CI, ihy2[iy]) );
  541. }
  542. // Third Off Diagonal
  543. if (iz!=nz-1) {
  544. triplets.push_back( Tr(ir , ic+(ny+1)*(nx) , -ahi_ihz2));
  545. triplets.push_back( Tr(ir+CI, ic+(ny+1)*(nx)+CI, ahi_ihz2));
  546. }
  547. ///////////////////////////////////////////////////////////////////////
  548. // imaginary part
  549. triplets.push_back( Tr(ir, ic+CI, omega*MU0*hsig) );
  550. ++ir;
  551. ++ic;
  552. }
  553. }
  554. }
  555. assert(ic == unx+uny);
  556. // LAPL{Az} + i omega mu sigmaz
  557. for (int iz=0; iz<nz+1; ++iz) {
  558. for (int iy=0; iy<ny; ++iy) {
  559. for (int ix=0; ix<nx; ++ix) {
  560. // Calculate 1/2 step values on staggered grid
  561. Real alo_ihx2(0);
  562. Real ahi_ihx2(0);
  563. if (ix == 0) {
  564. ahi_ihx2 = std::pow(.5*hx[ix+1] + .5*hx[ix], -2);
  565. alo_ihx2 = std::pow(.5*hx[ix+1] + .5*hx[ix], -2);
  566. } else if (ix == nx-1) {
  567. alo_ihx2 = std::pow(.5*hx[ix-1] + .5*hx[ix], -2);
  568. ahi_ihx2 = std::pow(.5*hx[ix-1] + .5*hx[ix], -2);
  569. } else {
  570. alo_ihx2 = std::pow(.5*hx[ix-1] + .5*hx[ix], -2);
  571. ahi_ihx2 = std::pow(.5*hx[ix+1] + .5*hx[ix], -2);
  572. }
  573. // Calculate 1/2 step values on staggered grid
  574. Real alo_ihy2(0);
  575. Real ahi_ihy2(0);
  576. if (iy == 0) {
  577. ahi_ihy2 = std::pow(.5*hy[iy+1] + .5*hy[iy], -2);
  578. alo_ihy2 = std::pow(.5*hy[iy+1] + .5*hy[iy], -2);
  579. } else if (iy == ny-1) {
  580. alo_ihy2 = std::pow(.5*hy[iy-1] + .5*hy[iy], -2);
  581. ahi_ihy2 = std::pow(.5*hy[iy-1] + .5*hy[iy], -2);
  582. } else {
  583. alo_ihy2 = std::pow(.5*hy[iy-1] + .5*hy[iy], -2);
  584. ahi_ihy2 = std::pow(.5*hy[iy+1] + .5*hy[iy], -2);
  585. }
  586. // Harmonically average sigmaz
  587. if (iz==0) {
  588. hsig = sigma[ix][iy][iz];
  589. } else if (iz == nz) {
  590. hsig = sigma[ix][iy][iz-1];
  591. } else {
  592. hsig = ((hz[iz]+hz[iz-1])/2.) *
  593. (1. / ( (hz[iz ] / (2.*sigma[ix][iy][iz ] + EPS)) +
  594. ( hz[iz-1] / (2.*sigma[ix][iy][iz-1] + EPS)) ) );
  595. }
  596. if (std::abs(hsig) < EPS2) hsig = 0;
  597. //hsig += Complex(0, omega*EPSILON0);
  598. Real Sum(0);
  599. if (iz == 0) {
  600. // Dirichlet bdry condition on A
  601. //Sum = alo_ihx2 + ahi_ihx2 + alo_ihy2 + ahi_ihy2 + 2.*ihz2[iz ];
  602. // Neumann bdry
  603. Sum = alo_ihx2 + ahi_ihx2 + alo_ihy2 + ahi_ihy2 + ihz2[iz ];
  604. } else if (iz == nz) {
  605. // Dirichlet bdry condition on A
  606. //Sum = alo_ihx2 + ahi_ihx2 + alo_ihy2 + ahi_ihy2 + 2.*ihz2[iz-1];
  607. // Neumann bdry
  608. Sum = alo_ihx2 + ahi_ihx2 + alo_ihy2 + ahi_ihy2 + ihz2[iz-1];
  609. } else {
  610. Sum = alo_ihx2 + ahi_ihx2 + alo_ihy2 + ahi_ihy2 + ihz2[iz] + ihz2[iz-1];
  611. }
  612. //////////////////////////////////////////////////////////////////
  613. // Diagonal Term
  614. triplets.push_back( Tr(ir , ic , Sum) );
  615. triplets.push_back( Tr(ir+CI, ic+CI, -Sum) );
  616. // First Off Diagonal
  617. if (ix!=nx-1) {
  618. triplets.push_back( Tr(ir , ic+1 , -ahi_ihx2) );
  619. triplets.push_back( Tr(ir+CI, ic+1+CI, ahi_ihx2) );
  620. }
  621. // Second Off Diagonal
  622. if (iy!=ny-1) {
  623. triplets.push_back( Tr(ir , ic+(nx) , -ahi_ihy2) );
  624. triplets.push_back( Tr(ir+CI, ic+(nx)+CI, ahi_ihy2) );
  625. }
  626. // Third Off Diagonal
  627. if (iz!=nz) {
  628. triplets.push_back( Tr(ir , ic+ny*(nx) , -ihz2[iz]) );
  629. triplets.push_back( Tr(ir+CI, ic+ny*(nx)+CI, ihz2[iz]) );
  630. }
  631. //////////////////////////////////////////////////////
  632. // imaginary part
  633. triplets.push_back( Tr(ir, ic+CI, omega*MU0*hsig) );
  634. ++ir;
  635. ++ic;
  636. }
  637. }
  638. }
  639. assert(ic == unx+uny+unz);
  640. CvecRe[iw].setFromTriplets(triplets.begin(), triplets.end()) ;
  641. CvecRe[iw].makeCompressed();
  642. return ;
  643. } // ----- end of method EMSchur3DBase::BuildCReal -----
  644. */
  645. //--------------------------------------------------------------------------------------
  646. // Class: EMSchur3DBase
  647. // Method: SetResFileName
  648. //--------------------------------------------------------------------------------------
  649. void EMSchur3DBase::SetResFileName ( const std::string& fname ) {
  650. ResFile = fname;
  651. return ;
  652. } // ----- end of method EMSchur3DBase::SetResFileName -----
  653. //--------------------------------------------------------------------------------------
  654. // Class: EMSchur3DBase
  655. // Method: BuildCPreconditioner
  656. //--------------------------------------------------------------------------------------
  657. // void EMSchur3DBase::BuildCPreconditioner ( const int& iw ) {
  658. // std::cout << "Building preconditioner for C_" << iw << std::endl;
  659. // CMvec[iw].setDroptol(1e-4); // 1e-4 previously (tested value)
  660. // CMvec[iw].setFillfactor(10);
  661. // Eigen::SparseMatrix<Complex> Csym;
  662. // Csym = Cvec[iw].selfadjointView<Eigen::Upper>();
  663. // CMvec[iw].compute( Csym );
  664. // return ;
  665. // } // ----- end of method EMSchur3DBase::BuildCPreconditioner -----
  666. //--------------------------------------------------------------------------------------
  667. // Class: EMSchur3DBase
  668. // Method: BuildD
  669. //--------------------------------------------------------------------------------------
  670. void EMSchur3DBase::BuildD ( ) {
  671. D.resize(uns, unx+uny+unz );
  672. //D.reserve(Eigen::VectorXi::Constant(unx+uny+unz, 6));
  673. std::vector<Tc> triplets;
  674. triplets.reserve( 2*unx + 2*uny + 2*unz );
  675. std::cout << "Building D...";
  676. int ic = 0;
  677. int ir = 0;
  678. ////////////////////////////////
  679. // Ax dx
  680. for (int iz=0; iz<nz; ++iz) {
  681. for (int iy=0; iy<ny; ++iy) {
  682. ++ic; // TRICKY
  683. for (int ix=0; ix<nx; ++ix) {
  684. //D.insert(ir, ic-1) = Complex(-ihx[ix], 0);
  685. //D.insert(ir, ic ) = Complex( ihx[ix], 0);
  686. triplets.push_back( Tc(ir, ic-1, Complex(-ihx[ix], 0)) );
  687. triplets.push_back( Tc(ir, ic , Complex( ihx[ix], 0)) );
  688. ++ir;
  689. ++ic;
  690. }
  691. }
  692. }
  693. assert (ic==unx);
  694. ///////////////////////////////
  695. // Ay dy
  696. ir = 0;
  697. for (int iz=0; iz<nz; ++iz) {
  698. for (int iy=0; iy<ny; ++iy) {
  699. for (int ix=0; ix<nx; ++ix) {
  700. //D.insert(ir, ic ) = Complex(-ihy[iy], 0);
  701. //D.insert(ir, ic+nx) = Complex( ihy[iy], 0);
  702. triplets.push_back( Tc(ir, ic , Complex(-ihy[iy], 0)) );
  703. triplets.push_back( Tc(ir, ic+nx, Complex( ihy[iy], 0)) );
  704. ++ir;
  705. ++ic;
  706. }
  707. }
  708. ic += nx;
  709. }
  710. assert (ic==unx+uny);
  711. ///////////////////////////////
  712. // Az dz
  713. ir = 0;
  714. for (int iz=0; iz<nz; ++iz) {
  715. for (int iy=0; iy<ny; ++iy) {
  716. for (int ix=0; ix<nx; ++ix) {
  717. //D.insert(ir, ic ) = Complex(-ihz[iz], 0);
  718. //D.insert(ir, ic+nx*ny) = Complex( ihz[iz], 0);
  719. triplets.push_back( Tc(ir, ic , Complex(-ihz[iz], 0)) );
  720. triplets.push_back( Tc(ir, ic+nx*ny, Complex( ihz[iz], 0)) );
  721. ++ir;
  722. ++ic;
  723. }
  724. }
  725. }
  726. assert (ic+nx*ny==unx+uny+unz);
  727. D.setFromTriplets(triplets.begin(), triplets.end());
  728. D.makeCompressed();
  729. std::cout << "Done!" << std::endl;
  730. return ;
  731. } // ----- end of method EMSchur3DBase::BuildD -----
  732. //--------------------------------------------------------------------------------------
  733. // Class: EMSchur3DBase
  734. // Method: SetRectilinearGrid
  735. //--------------------------------------------------------------------------------------
  736. void EMSchur3DBase::SetRectilinearGrid ( std::shared_ptr<RectilinearGrid> GridPtr ) {
  737. Grid = GridPtr;
  738. // Do some convienience stuff, we call these so often it's nice to not have to use
  739. // the accessors every time. This is the only place these are set.
  740. nx = Grid->GetNx();
  741. ny = Grid->GetNy();
  742. nz = Grid->GetNz();
  743. hx = Grid->GetDx();
  744. hy = Grid->GetDy();
  745. hz = Grid->GetDz();
  746. unx = (nx+1)*ny*nz; // Number of Vectors x component
  747. uny = nx*(ny+1)*nz; // Number of Vectors y component
  748. unz = nx*ny*(nz+1); // Number of Vectors z component
  749. uns = nx*ny*nz; // On dual grid, number of scalars
  750. ihx = 1. / hx.array();
  751. ihy = 1. / hy.array();
  752. ihz = 1. / hz.array();
  753. ihx2 = 1. / (hx.array() * hx.array());
  754. ihy2 = 1. / (hy.array() * hy.array());
  755. ihz2 = 1. / (hz.array() * hz.array());
  756. return ;
  757. } // ----- end of method EMSchur3DBase::SetRectilinearGrid -----
  758. //--------------------------------------------------------------------------------------
  759. // Class: EMSchur3DBase
  760. // Method: SetRectilinearGrid
  761. //--------------------------------------------------------------------------------------
  762. void EMSchur3DBase::SetRectilinearGrid ( vtkRectilinearGrid* rGridPtr ) {
  763. // Convert to Lemma here
  764. } // ----- end of method EMSchur3DBase::SetRectilinearGrid -----
  765. //--------------------------------------------------------------------------------------
  766. // Class: EMSchur3DBase
  767. // Method: Setup
  768. //--------------------------------------------------------------------------------------
  769. void EMSchur3DBase::Setup ( ) {
  770. ////////////////////////////////////
  771. // First set up grid stuff
  772. if (Cvec) delete [] Cvec;
  773. //if (CSolver) delete [] CSolver;
  774. Cvec = new Eigen::SparseMatrix<Complex, Eigen::ColMajor> [Omegas.size()];
  775. //#ifdef LEMMAUSEOMP
  776. //#pragma omp parallel for schedule(static, 1)
  777. //#endif
  778. for (int iw=0; iw<Omegas.size(); ++iw) {
  779. std::cout << "Build C " << std::endl;
  780. BuildC(sigma, sigma, sigma, iw);
  781. // std::cout << "Build Re(C) " << std::endl;
  782. // BuildCReal(sigma, sigma, sigma, iw);
  783. }
  784. BuildCDirectSolver( ); // templated specialisation
  785. BuildD();
  786. return ;
  787. } // ----- end of method EMSchur3DBase::Setup -----
  788. //--------------------------------------------------------------------------------------
  789. // Class: EMSchur3DBase
  790. // Method: SetLayeredEarthEM
  791. //--------------------------------------------------------------------------------------
  792. void EMSchur3DBase::SetLayeredEarthEM ( std::shared_ptr<LayeredEarthEM> LayModelin ) {
  793. LayModel = LayModelin;
  794. return ;
  795. } // ----- end of method EMSchur3DBase::SetLayeredEarthEM -----
  796. //--------------------------------------------------------------------------------------
  797. // Class: EMSchur3DBase
  798. // Method: PrimaryField
  799. //--------------------------------------------------------------------------------------
  800. void EMSchur3DBase::PrimaryField ( std::shared_ptr<DipoleSource> source, std::shared_ptr<FieldPoints> dpoint ) {
  801. auto EmEarth = Lemma::EMEarth1D::NewSP();
  802. EmEarth->AttachDipoleSource(source);
  803. EmEarth->AttachLayeredEarthEM(LayModel);
  804. EmEarth->AttachFieldPoints(dpoint);
  805. EmEarth->SetFieldsToCalculate(Lemma::E);
  806. EmEarth->SetHankelTransformMethod(ANDERSON801);
  807. std::cout << "Primary Field Calculation" << std::endl;
  808. EmEarth->MakeCalc3();
  809. return ;
  810. } // ----- end of method EMSchur3DBase::PrimaryField -----
  811. //--------------------------------------------------------------------------------------
  812. // Class: EMSchur3DBase
  813. // Method: FillPoints
  814. //--------------------------------------------------------------------------------------
  815. void EMSchur3DBase::FillPoints ( std::shared_ptr<FieldPoints> dpoint) {
  816. dpoint->SetNumberOfPoints(unx + uny + unz);
  817. int ip(0);
  818. Real tx, ty, tz;
  819. Real ox = Grid->GetOx();
  820. Real oy = Grid->GetOy();
  821. Real oz = Grid->GetOz();
  822. tz = oz;
  823. for (int iz=0; iz<nz; ++iz) {
  824. ty = oy;
  825. for (int iy=0; iy<ny; ++iy) {
  826. tx = ox - hx(0)/2.;
  827. for (int ix=0; ix<nx+1; ++ix) {
  828. dpoint->SetLocation(ip, tx, ty, tz);
  829. if (ix < nx) tx += hx(ix);
  830. ++ip;
  831. }
  832. if (iy<ny-1) ty += .5*hy(iy) + .5*hy(iy+1);
  833. }
  834. if (iz<nz-1) tz += .5*hz(iz) + .5*hz(iz+1);
  835. }
  836. tz = oz;
  837. for (int iz=0; iz<nz; ++iz) {
  838. ty = oy - hy(0)/2.;
  839. for (int iy=0; iy<ny+1; ++iy) {
  840. tx = ox;
  841. for (int ix=0; ix<nx; ++ix) {
  842. dpoint->SetLocation(ip, tx, ty, tz);
  843. if (ix<nx-1) tx += .5*hx(ix) + .5*hx(ix+1);
  844. ++ip;
  845. }
  846. if (iy < ny) ty += hy(iy);
  847. }
  848. if (iz<nz-1) tz += .5*hz(iz) + .5*hz(iz+1);
  849. }
  850. tz = oz - hz(0)/2.;
  851. for (int iz=0; iz<nz+1; ++iz) {
  852. ty = oy;
  853. for (int iy=0; iy<ny; ++iy) {
  854. tx = ox;
  855. for (int ix=0; ix<nx; ++ix) {
  856. dpoint->SetLocation(ip, tx, ty, tz);
  857. if (ix<nx-1) tx += .5*hx(ix) + .5*hx(ix+1);
  858. ++ip;
  859. }
  860. if (iy<ny-1) ty += .5*hy(iy) + .5*hy(iy+1);
  861. }
  862. if (iz < nz) tz += hz(iz);
  863. }
  864. return ;
  865. } // ----- end of method EMSchur3DBase::FillPoints -----
  866. //--------------------------------------------------------------------------------------
  867. // Class: EMSchur3DBase
  868. // Method: FillSourceTerms
  869. //--------------------------------------------------------------------------------------
  870. void EMSchur3DBase::FillSourceTerms ( Eigen::Ref<VectorXcr> ms,
  871. Eigen::Ref<VectorXcr> Se, Eigen::Ref<VectorXcr> E0,
  872. std::shared_ptr<FieldPoints> dpoint, const Real& omega ) {
  873. //Complex hsig(0);
  874. //Complex hsigp(0);
  875. Real hsig(0);
  876. Real hsigp(0);
  877. Real EPS(0); //1e-24);
  878. Real EPS2(1e-18);
  879. int iv = 0;
  880. for (int iz=0; iz<nz; ++iz) {
  881. for (int iy=0; iy<ny; ++iy) {
  882. for (int ix=0; ix<nx+1; ++ix) {
  883. // Harmonically average sigma
  884. if (ix == 0) {
  885. hsig = sigma[ix][iy][iz];
  886. hsigp = sigmap[ix][iy][iz];
  887. } else if (ix == nx) {
  888. hsig = sigma[ix-1][iy][iz];
  889. hsigp = sigmap[ix-1][iy][iz];
  890. } else {
  891. hsig = ((hx[ix]+hx[ix-1])/2.) *
  892. (1. / ( (hx[ix ] / (2.*sigma[ix ][iy][iz] + EPS)) +
  893. ( hx[ix-1] / (2.*sigma[ix-1][iy][iz] + EPS)) ) );
  894. hsigp = ((hx[ix]+hx[ix-1])/2.) *
  895. (1. / ( (hx[ix ] / (2.*sigmap[ix ][iy][iz] + EPS)) +
  896. ( hx[ix-1] / (2.*sigmap[ix-1][iy][iz] + EPS)) ) );
  897. }
  898. if (std::abs(hsig) < EPS2) hsig = 0;
  899. if (std::abs(hsigp) < EPS2) hsigp = 0;
  900. //hsig += Complex(0, omega*EPSILON0);
  901. //ioms(iv) = Complex(0, omega*MU0*hsig);
  902. //Se(iv) = ioms(iv)*dpoint->GetEfield(0,iv)(0);
  903. ms(iv) = MU0*hsig;
  904. Se(iv) = MU0*hsigp*
  905. dpoint->GetEfield(0,iv)(0);
  906. E0(iv) = dpoint->GetEfield(0,iv)(0);
  907. ++iv;
  908. }
  909. }
  910. }
  911. for (int iz=0; iz<nz; ++iz) {
  912. for (int iy=0; iy<ny+1; ++iy) {
  913. for (int ix=0; ix<nx; ++ix) {
  914. // Harmonically average sigma
  915. if (iy == 0) {
  916. hsig = sigma[ix][iy][iz];
  917. hsigp = sigmap[ix][iy][iz];
  918. } else if (iy == ny) {
  919. hsig = sigma[ix][iy-1][iz];
  920. hsigp = sigmap[ix][iy-1][iz];
  921. } else {
  922. hsig = ((hy[iy]+hy[iy-1])/2.) *
  923. (1. / ( (hy[iy ] / (2.*sigma[ix][iy ][iz] + EPS)) +
  924. ( hy[iy-1] / (2.*sigma[ix][iy-1][iz] + EPS)) ) );
  925. hsigp = ((hy[iy]+hy[iy-1])/2.) *
  926. (1. / ( (hy[iy ] / (2.*sigmap[ix][iy ][iz] + EPS)) +
  927. ( hy[iy-1] / (2.*sigmap[ix][iy-1][iz] + EPS)) ) );
  928. }
  929. if (std::abs(hsig) < EPS2) hsig = 0;
  930. if (std::abs(hsigp) < EPS2) hsigp = 0;
  931. //hsig += Complex(0, omega*EPSILON0);
  932. //ioms(iv) = Complex(0, omega*MU0*hsig);
  933. //Se(iv) = ioms(iv) * dpoint->GetEfield(0,iv)(1);
  934. ms(iv) = MU0*hsig;
  935. Se(iv) = MU0*hsigp* // TODO was hsigp !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  936. dpoint->GetEfield(0,iv)(1);
  937. E0(iv) = dpoint->GetEfield(0,iv)(1);
  938. ++iv;
  939. }
  940. }
  941. }
  942. for (int iz=0; iz<nz+1; ++iz) {
  943. for (int iy=0; iy<ny; ++iy) {
  944. for (int ix=0; ix<nx; ++ix) {
  945. // Harmonically average sigma
  946. if (iz==0) {
  947. hsig = sigma[ix][iy][iz];
  948. hsigp = sigma[ix][iy][iz];
  949. } else if (iz == nz) {
  950. hsig = sigma[ix][iy][nz-1];
  951. hsigp = sigmap[ix][iy][nz-1];
  952. } else {
  953. hsig = ((hz[iz]+hz[iz-1])/2.) *
  954. (1. / ( (hz[iz ] / (2.*sigma[ix][iy][iz ] + EPS)) +
  955. ( hz[iz-1] / (2.*sigma[ix][iy][iz-1] + EPS)) ) );
  956. hsigp = ((hz[iz]+hz[iz-1])/2.) *
  957. (1. / ( (hz[iz ] / (2.*sigmap[ix][iy][iz ] + EPS)) +
  958. ( hz[iz-1] / (2.*sigmap[ix][iy][iz-1] + EPS)) ) );
  959. }
  960. if (std::abs(hsig) < EPS2) hsig = 0;
  961. if (std::abs(hsigp) < EPS2) hsigp = 0;
  962. //hsig += Complex(0, omega*EPSILON0);
  963. //ioms(iv) = Complex(0, omega*MU0*hsig);
  964. //Se(iv) = ioms(iv) * dpoint->GetEfield(0, iv)(2);
  965. ms(iv) = MU0*hsig;
  966. Se(iv) = MU0*hsig
  967. *dpoint->GetEfield(0, iv)(2);
  968. E0(iv) = dpoint->GetEfield(0, iv)(2);
  969. ++iv;
  970. }
  971. }
  972. }
  973. return ;
  974. } // ----- end of method EMSchur3DBase::FillSourceTerms -----
  975. VectorXcr EMSchur3DBase::StaggeredGridCurl( Eigen::Ref<VectorXcr> A ) {
  976. VectorXcr B = VectorXcr::Zero(3*nx*ny*nz);
  977. VectorXcr dzdy = VectorXcr::Zero( nx*ny*nz);
  978. VectorXcr dydz = VectorXcr::Zero( nx*ny*nz);
  979. VectorXcr dxdz = VectorXcr::Zero( nx*ny*nz);
  980. VectorXcr dzdx = VectorXcr::Zero( nx*ny*nz);
  981. VectorXcr dydx = VectorXcr::Zero( nx*ny*nz);
  982. VectorXcr dxdy = VectorXcr::Zero( nx*ny*nz);
  983. /*
  984. Curl_x = dzdy - dydz
  985. Curl_y = dxdz - dzdx
  986. Curl_z = dydx - dxdy
  987. */
  988. ////////////////////////////////
  989. // X component
  990. int ii = 0;
  991. int iix = 0;
  992. for (int iz=0; iz<nz; ++iz) {
  993. for (int iy=0; iy<ny; ++iy) {
  994. for (int ix=0; ix<nx; ++ix) {
  995. dxdz(iix) = ( A(ii+(nx+1)*ny) - A(ii) ) / hz[iz];
  996. dxdy(iix) = ( A(ii+(nx+1) ) - A(ii) ) / hy[iy];
  997. ++ii;
  998. ++iix;
  999. }
  1000. // Jump around boundary
  1001. ++ii;
  1002. }
  1003. }
  1004. ///////////////////////////////
  1005. // Y component
  1006. int iiy = 0;
  1007. for (int iz=0; iz<nz; ++iz) {
  1008. for (int iy=0; iy<ny; ++iy) {
  1009. for (int ix=0; ix<nx; ++ix) {
  1010. dydz(iiy) = ( A(ii+(nx)*(ny+1)) - A(ii) ) / hz[iz];
  1011. dydx(iiy) = ( A(ii+1 ) - A(ii) ) / hx[ix];
  1012. ++iiy;
  1013. ++ii;
  1014. }
  1015. }
  1016. // Jump around boundary
  1017. ii += nx;
  1018. }
  1019. ////////////////////////////
  1020. // Z component
  1021. int iiz = 0;
  1022. for (int iz=0; iz<nz; ++iz) {
  1023. for (int iy=0; iy<ny; ++iy) {
  1024. for (int ix=0; ix<nx; ++ix) {
  1025. //Az(iiz) = 0.5 * ( A(ii) + A(ii+nx*ny));
  1026. dzdy(iiz) = ( A(ii+nx ) - A(ii) ) / hy[iy];
  1027. dzdx(iiz) = ( A(ii+1 ) - A(ii) ) / hx[ix];
  1028. ++iiz;
  1029. ++ii;
  1030. }
  1031. //++iiz;
  1032. //ii += nx*ny;
  1033. }
  1034. }
  1035. B.segment( 0, nx*ny*nz) = dzdy - dydz;
  1036. B.segment( nx*ny*nz, nx*ny*nz) = dxdz - dzdx;
  1037. B.segment(2*nx*ny*nz, nx*ny*nz) = dydx - dxdy;
  1038. return B;
  1039. }
  1040. //--------------------------------------------------------------------------------------
  1041. // Class: EMSchur3DBase
  1042. // Method: WriteVTKResults
  1043. //--------------------------------------------------------------------------------------
  1044. void EMSchur3DBase::WriteVTKResults ( const std::string& fname, Eigen::Ref<VectorXcr> A,
  1045. Eigen::Ref<VectorXcr> Se, Eigen::Ref<VectorXcr> E0, Eigen::Ref<VectorXcr> E,
  1046. Eigen::Ref<VectorXcr> Phi, Eigen::Ref<VectorXcr> ADiv, Eigen::Ref<VectorXcr> ADiv2,
  1047. Eigen::Ref<VectorXcr> B ) {
  1048. auto VTKGridExporter = RectilinearGridVTKExporter::NewSP();
  1049. VTKGridExporter->SetGrid(Grid);
  1050. // VTK arrays to fill
  1051. // scalars
  1052. vtkDoubleArray *sigmaArray = vtkDoubleArray::New();
  1053. vtkDoubleArray *sigmapArray = vtkDoubleArray::New();
  1054. vtkDoubleArray *phiRealArray = vtkDoubleArray::New();
  1055. vtkDoubleArray *phiImagArray = vtkDoubleArray::New();
  1056. vtkDoubleArray *ADivRealArray = vtkDoubleArray::New();
  1057. vtkDoubleArray *ADivImagArray = vtkDoubleArray::New();
  1058. vtkDoubleArray *ADiv2RealArray = vtkDoubleArray::New();
  1059. vtkDoubleArray *ADiv2ImagArray = vtkDoubleArray::New();
  1060. // vectors
  1061. vtkDoubleArray *AReal = vtkDoubleArray::New();
  1062. vtkDoubleArray *AImag = vtkDoubleArray::New();
  1063. vtkDoubleArray *BReal = vtkDoubleArray::New();
  1064. vtkDoubleArray *BImag = vtkDoubleArray::New();
  1065. vtkDoubleArray *SeReal = vtkDoubleArray::New();
  1066. vtkDoubleArray *SeImag = vtkDoubleArray::New();
  1067. vtkDoubleArray *E0Real = vtkDoubleArray::New();
  1068. vtkDoubleArray *E0Imag = vtkDoubleArray::New();
  1069. vtkDoubleArray *EReal = vtkDoubleArray::New();
  1070. vtkDoubleArray *EImag = vtkDoubleArray::New();
  1071. AReal->SetNumberOfComponents(3);
  1072. AImag->SetNumberOfComponents(3);
  1073. BReal->SetNumberOfComponents(3);
  1074. BImag->SetNumberOfComponents(3);
  1075. SeReal->SetNumberOfComponents(3);
  1076. SeImag->SetNumberOfComponents(3);
  1077. E0Real->SetNumberOfComponents(3);
  1078. E0Imag->SetNumberOfComponents(3);
  1079. EReal->SetNumberOfComponents(3);
  1080. EImag->SetNumberOfComponents(3);
  1081. // Interpolate
  1082. Eigen::VectorXcd Ax(nx*ny*nz); // A
  1083. Eigen::VectorXcd Ay(nx*ny*nz);
  1084. Eigen::VectorXcd Az(nx*ny*nz);
  1085. Eigen::VectorXcd Bx(nx*ny*nz); // B
  1086. Eigen::VectorXcd By(nx*ny*nz);
  1087. Eigen::VectorXcd Bz(nx*ny*nz);
  1088. Eigen::VectorXcd Sex(nx*ny*nz); // Se
  1089. Eigen::VectorXcd Sey(nx*ny*nz);
  1090. Eigen::VectorXcd Sez(nx*ny*nz);
  1091. Eigen::VectorXcd E0x(nx*ny*nz); // E0
  1092. Eigen::VectorXcd E0y(nx*ny*nz);
  1093. Eigen::VectorXcd E0z(nx*ny*nz);
  1094. Eigen::VectorXcd Ex(nx*ny*nz); // E
  1095. Eigen::VectorXcd Ey(nx*ny*nz);
  1096. Eigen::VectorXcd Ez(nx*ny*nz);
  1097. ////////////////////////////////
  1098. // X component
  1099. int ii = 0;
  1100. int iix = 0;
  1101. int ic = 0;
  1102. for (int iz=0; iz<nz; ++iz) {
  1103. for (int iy=0; iy<ny; ++iy) {
  1104. for (int ix=0; ix<nx; ++ix) {
  1105. Ax(iix) = 0.5 * ( A(ii) + A(ii+1));
  1106. Sex(iix) = 0.5 * (Se(ii) + Se(ii+1));
  1107. E0x(iix) = 0.5 * (E0(ii) + E0(ii+1));
  1108. Ex(iix) = 0.5 * ( E(ii) + E(ii+1));
  1109. Bx(iix) = B(ic);
  1110. ++ii;
  1111. ++iix;
  1112. ++ic;
  1113. }
  1114. // Jump around boundary
  1115. ++ii;
  1116. }
  1117. }
  1118. ///////////////////////////////
  1119. // Y component
  1120. int iiy = 0;
  1121. for (int iz=0; iz<nz; ++iz) {
  1122. for (int iy=0; iy<ny; ++iy) {
  1123. for (int ix=0; ix<nx; ++ix) {
  1124. Ay(iiy) = 0.5 * ( A(ii) + A(ii+nx));
  1125. Sey(iiy) = 0.5 * (Se(ii) + Se(ii+nx));
  1126. E0y(iiy) = 0.5 * (E0(ii) + E0(ii+nx));
  1127. Ey(iiy) = 0.5 * ( E(ii) + E(ii+nx));
  1128. By(iiy) = B(ic);
  1129. ++iiy;
  1130. ++ii;
  1131. ++ic;
  1132. }
  1133. }
  1134. // Jump around boundary
  1135. ii += nx;
  1136. }
  1137. ////////////////////////////
  1138. // Z component
  1139. int iiz = 0;
  1140. for (int iz=0; iz<nz; ++iz) {
  1141. for (int iy=0; iy<ny; ++iy) {
  1142. for (int ix=0; ix<nx; ++ix) {
  1143. Az(iiz) = 0.5 * ( A(ii) + A(ii+nx*ny));
  1144. Sez(iiz) = 0.5 * (Se(ii) + Se(ii+nx*ny));
  1145. E0z(iiz) = 0.5 * (E0(ii) + E0(ii+nx*ny));
  1146. Ez(iiz) = 0.5 * ( E(ii) + E(ii+nx*ny));
  1147. Bz(iiz) = B(ic);
  1148. ++iiz;
  1149. ++ii;
  1150. ++ic;
  1151. }
  1152. //++iiz;
  1153. //ii += nx*ny;
  1154. }
  1155. }
  1156. //VTKGridExporter->GetVTKGrid();
  1157. int i = 0;
  1158. for (int iz=0; iz<nz; ++iz) {
  1159. for (int iy=0; iy<ny; ++iy) {
  1160. for (int ix=0; ix<nx; ++ix) {
  1161. // scalars
  1162. //sigmaArray->InsertTuple1(i, std::real(sigma[ix][iy][iz]) );
  1163. //sigmapArray->InsertTuple1(i, std::real(sigmap[ix][iy][iz]) );
  1164. sigmaArray->InsertTuple1(i, sigma[ix][iy][iz] );
  1165. sigmapArray->InsertTuple1(i, sigmap[ix][iy][iz] );
  1166. phiRealArray->InsertTuple1(i, std::real(Phi(i)));
  1167. phiImagArray->InsertTuple1(i, std::imag(Phi(i)));
  1168. ADivRealArray->InsertTuple1(i, std::real(ADiv(i)));
  1169. ADivImagArray->InsertTuple1(i, std::imag(ADiv(i)));
  1170. ADiv2RealArray->InsertTuple1(i, std::real(ADiv2(i)));
  1171. ADiv2ImagArray->InsertTuple1(i, std::imag(ADiv2(i)));
  1172. // vectors
  1173. AReal->InsertTuple3(i, std::real(Ax(i)), std::real(Ay(i)), std::real(Az(i)));
  1174. AImag->InsertTuple3(i, std::imag(Ax(i)), std::imag(Ay(i)), std::imag(Az(i)));
  1175. BReal->InsertTuple3(i, std::real(Bx(i)), std::real(By(i)), std::real(Bz(i)));
  1176. BImag->InsertTuple3(i, std::imag(Bx(i)), std::imag(By(i)), std::imag(Bz(i)));
  1177. SeReal->InsertTuple3(i, std::real(Sex(i)), std::real(Sey(i)), std::real(Sez(i)));
  1178. SeImag->InsertTuple3(i, std::imag(Sex(i)), std::imag(Sey(i)), std::imag(Sez(i)));
  1179. E0Real->InsertTuple3(i, std::real(E0x(i)), std::real(E0y(i)), std::real(E0z(i)));
  1180. E0Imag->InsertTuple3(i, std::imag(E0x(i)), std::imag(E0y(i)), std::imag(E0z(i)));
  1181. EReal-> InsertTuple3(i, std::real(Ex(i)), std::real(Ey(i)), std::real(Ez(i)));
  1182. EImag-> InsertTuple3(i, std::imag(Ex(i)), std::imag(Ey(i)), std::imag(Ez(i)));
  1183. ++i;
  1184. }
  1185. }
  1186. }
  1187. AReal->SetName("A_real");
  1188. AImag->SetName("A_imag");
  1189. BReal->SetName("B_real");
  1190. BImag->SetName("B_imag");
  1191. SeReal->SetName("Se_real");
  1192. SeImag->SetName("Se_imag");
  1193. E0Real->SetName("E0_real");
  1194. E0Imag->SetName("E0_imag");
  1195. EReal->SetName("E_real");
  1196. EImag->SetName("E_imag");
  1197. phiRealArray->SetName("phi Real");
  1198. phiImagArray->SetName("phi imag");
  1199. ADivRealArray->SetName("ini div A real");
  1200. ADivImagArray->SetName("ini div A imag");
  1201. ADiv2RealArray->SetName("div A real");
  1202. ADiv2ImagArray->SetName("div A imag");
  1203. // Add scalar fields
  1204. sigmaArray->SetName("sigma");
  1205. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(sigmaArray);
  1206. sigmapArray->SetName("sigma prime");
  1207. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(sigmapArray);
  1208. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(phiRealArray);
  1209. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(phiImagArray);
  1210. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(ADivRealArray);
  1211. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(ADivImagArray);
  1212. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(ADiv2RealArray);
  1213. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(ADiv2ImagArray);
  1214. // Add Vector Arrays
  1215. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(AReal);
  1216. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(AImag);
  1217. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(SeReal);
  1218. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(SeImag);
  1219. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(E0Real);
  1220. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(E0Imag);
  1221. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(EReal);
  1222. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(EImag);
  1223. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(BReal);
  1224. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(BImag);
  1225. // Write
  1226. vtkXMLRectilinearGridWriter *gridWrite = vtkXMLRectilinearGridWriter::New();
  1227. gridWrite->SetInputData( VTKGridExporter->GetVTKGrid() );
  1228. gridWrite->SetFileName( (fname + std::string(".vtr")).c_str());
  1229. gridWrite->Update();
  1230. gridWrite->Write();
  1231. sigmaArray->Delete();
  1232. gridWrite->Delete();
  1233. phiRealArray->Delete();
  1234. phiImagArray->Delete();
  1235. ADivRealArray->Delete();
  1236. ADivImagArray->Delete();
  1237. ADiv2RealArray->Delete();
  1238. ADiv2ImagArray->Delete();
  1239. // TODO delete vectors too!
  1240. AReal->Delete();
  1241. AImag->Delete();
  1242. BReal->Delete();
  1243. BImag->Delete();
  1244. SeReal->Delete();
  1245. SeImag->Delete();
  1246. E0Real->Delete();
  1247. E0Imag->Delete();
  1248. EReal->Delete();
  1249. EImag->Delete();
  1250. return ;
  1251. } // ----- end of method EMSchur3DBase::WriteVTKResults -----
  1252. //--------------------------------------------------------------------------------------
  1253. // Class: EMSchur3DBase
  1254. // Method: LoadMeshToolsConductivityModel
  1255. //--------------------------------------------------------------------------------------
  1256. void EMSchur3DBase::LoadMeshToolsConductivityModel ( const std::string& fname ) {
  1257. if (sigma) Delete3DScalar(sigma);
  1258. //Allocate3DScalar(sigma, Complex(0));
  1259. Allocate3DScalar(sigma, 0.);
  1260. // TODO, there are smarter ways to do this, since the bacground is 1D. But for now, we just have
  1261. // them the same. Eases logic in harmonic averaging. And it's just in one thread, so not so
  1262. // hateful
  1263. if (!LayModel) {
  1264. std::cerr << "In EMSchur3DBase::LoadMeshToolsConductivityModel, you need to set your bacground 1D model\n";
  1265. exit(EXIT_FAILURE);
  1266. }
  1267. if (sigmap) Delete3DScalar(sigmap);
  1268. //Allocate3DScalar(sigmap, Complex(0));
  1269. Allocate3DScalar(sigmap, 0.);
  1270. auto Parser = ASCIIParser::NewSP();
  1271. Parser->Open(fname);
  1272. std::vector<Real> sigmod = Parser->ReadReals( -1 );
  1273. if ( static_cast<int>(sigmod.size()) != nx*ny*nz) {
  1274. std::cerr << "In EMSchur3DBase::LoadMeshToolsConductivityModel model sizes are not in agreement\n";
  1275. exit(EXIT_FAILURE);
  1276. }
  1277. Parser->Close();
  1278. int ic(0);
  1279. for (int ix=0; ix<nx; ++ix) {
  1280. for (int iy=0; iy<ny; ++iy) {
  1281. Real pz = Grid->GetOz();
  1282. for (int iz=0; iz<nz; ++iz) {
  1283. if ( sigmod[ic] != sigmod[ic] ) {
  1284. std::cout << "sigmod problems " << ic << std::endl;
  1285. }
  1286. sigma[ix][iy][iz] = sigmod[ic];
  1287. //sigmap[ix][iy][iz] = sigmod[ic] - LayModel->GetLayerConductivity(LayModel->GetLayerAtThisDepth(pz));
  1288. sigmap[ix][iy][iz] = sigmod[ic] - LayModel->GetLayerConductivity(LayModel->GetLayerAtThisDepth(pz)).real() ;
  1289. ++ic;
  1290. pz += hz[iz];
  1291. }
  1292. }
  1293. }
  1294. // Marker and cell
  1295. MAC = VectorXi::Zero(nx*ny*nz);
  1296. idx.clear();
  1297. ic = 0;
  1298. for (int iz=0; iz<nz; ++iz) {
  1299. for (int iy=0; iy<ny; ++iy) {
  1300. for (int ix=0; ix<nx; ++ix) {
  1301. if (sigma[ix][iy][iz] > 0.) {
  1302. MAC(ic) = 1;
  1303. idx.push_back(ic);
  1304. }
  1305. ++ic;
  1306. }
  1307. }
  1308. }
  1309. return ;
  1310. }
  1311. //--------------------------------------------------------------------------------------
  1312. // Class: EMSchur3DBase
  1313. // Method: SetAEMSurvey
  1314. //--------------------------------------------------------------------------------------
  1315. void EMSchur3DBase::SetAEMSurvey ( std::shared_ptr<AEMSurvey> SurveyIn ) {
  1316. Survey = SurveyIn;
  1317. // determine how many frequencies we are dealing with
  1318. Omegas = 2.*PI*Survey->GetFrequencies();
  1319. return ;
  1320. } // ----- end of method EMSchur3DBase::SetAEMSurvey -----// ----- end of method EMSchur3DBase::LoadMeshToolsConductivityModel -----
  1321. //--------------------------------------------------------------------------------------
  1322. // Class: EMSchur3DBase
  1323. // Method: Solve
  1324. //--------------------------------------------------------------------------------------
  1325. void EMSchur3DBase::Solve ( ) {
  1326. // TODO, these should throw exceptions for a GUI to catch
  1327. if (!LayModel or !Survey or !Grid or !sigma) {
  1328. std::cerr << "You need to set something. EMSChur3D::Solve()";
  1329. exit(EXIT_FAILURE);
  1330. }
  1331. //BuildD(); // strangely necessary, bug, track down.
  1332. Setup();
  1333. std::cout << "Entering parallel solve" << std::endl;
  1334. //#ifdef LEMMAUSEOMP
  1335. //#pragma omp parallel for schedule(static,1)
  1336. //#endif
  1337. for (int isource=0; isource<Survey->GetNumberOfSources(); ++isource) {
  1338. SolveSource(Survey->GetSource(isource), isource);
  1339. }
  1340. return ;
  1341. } // ----- end of method EMSchur3DBase::Solve -----
  1342. } // ----- end of Lemma name -----