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

EMSchur3DBase.cpp 61KB

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