3D EM based on Schur decomposition
Вы не можете выбрать более 25 тем Темы должны начинаться с буквы или цифры, могут содержать дефисы(-) и должны содержать не более 35 символов.

EMSchur3DBase.cpp 59KB

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