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

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