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 58KB

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