3D EM based on Schur decomposition
You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

EMSchur3DBase.cpp 59KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475147614771478147914801481148214831484148514861487148814891490149114921493149414951496149714981499150015011502150315041505150615071508150915101511151215131514
  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";
  31. return stream;
  32. }
  33. // ==================== LIFECYCLE =======================
  34. //--------------------------------------------------------------------------------------
  35. // Class: EMSchur3DBase
  36. // Method: EMSchur3DBase
  37. // Description: constructor (protected)
  38. //--------------------------------------------------------------------------------------
  39. // std::shared_ptr<EMSchur3DBase> EMSchur3DBase::NewSP() {
  40. // return std::make_shared<EMSchur3DBase>( ctor_key() );
  41. // }
  42. //--------------------------------------------------------------------------------------
  43. // Class: EMSchur3DBase
  44. // Method: EMSchur3DBase
  45. // Description: constructor (protected)
  46. //--------------------------------------------------------------------------------------
  47. EMSchur3DBase::EMSchur3DBase ( const ctor_key& key ) : LemmaObject( key ),
  48. Grid(nullptr),
  49. //Survey(nullptr),
  50. LayModel(nullptr), Cvec(nullptr),
  51. ResFile("source"), sigma(nullptr), sigmap(nullptr)
  52. {
  53. } // ----- end of method EMSchur3DBase::EMSchur3DBase (constructor) -----
  54. //--------------------------------------------------------------------------------------
  55. // Class: EMSchur3DBase
  56. // Method: EMSchur3DBase
  57. // Description: constructor (protected)
  58. //--------------------------------------------------------------------------------------
  59. EMSchur3DBase::EMSchur3DBase ( const YAML::Node& node, const ctor_key& key ) : LemmaObject( key ),
  60. Grid(nullptr),
  61. //Survey(nullptr),
  62. LayModel(nullptr), Cvec(nullptr),
  63. ResFile("source"), sigma(nullptr), sigmap(nullptr)
  64. {
  65. } // ----- end of method EMSchur3DBase::~EMSchur3DBase (destructor) -----
  66. //--------------------------------------------------------------------------------------
  67. // Class: EMSchur3DBase
  68. // Method: ~EMSchur3DBase
  69. // Description: destructor (protected)
  70. //--------------------------------------------------------------------------------------
  71. EMSchur3DBase::~EMSchur3DBase ( ) {
  72. if (sigma) Delete3DScalar(sigma);
  73. if (sigmap) Delete3DScalar(sigmap);
  74. if (Cvec) delete [] Cvec;
  75. //if (CvecRe) delete [] CvecRe;
  76. } // ----- end of method EMSchur3DBase::~EMSchur3DBase (destructor) -----
  77. //--------------------------------------------------------------------------------------
  78. // Class: EMSchur3D
  79. // Method: Serialize
  80. //--------------------------------------------------------------------------------------
  81. YAML::Node EMSchur3DBase::Serialize ( ) const {
  82. YAML::Node node = LemmaObject::Serialize();
  83. node.SetTag( this->GetName() );
  84. node["EMSchur3D_VERSION"] = EMSCHUR3D_VERSION;
  85. return node;
  86. } // ----- end of method EMSchur3D::Serialize -----
  87. //--------------------------------------------------------------------------------------
  88. // Class: EMSchur3DBase
  89. // Method: BuildC
  90. //--------------------------------------------------------------------------------------
  91. void EMSchur3DBase::BuildC ( Real*** sigmax, Real*** sigmay, Real*** sigmaz, const int& iw) {
  92. Cvec[iw].resize( unx+uny+unz , unx+uny+unz );
  93. 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: SetRectilinearGrid
  758. //--------------------------------------------------------------------------------------
  759. void EMSchur3DBase::SetRectilinearGrid ( vtkRectilinearGrid* rGridPtr ) {
  760. // Convert to Lemma here
  761. } // ----- end of method EMSchur3DBase::SetRectilinearGrid -----
  762. //--------------------------------------------------------------------------------------
  763. // Class: EMSchur3DBase
  764. // Method: Setup
  765. //--------------------------------------------------------------------------------------
  766. void EMSchur3DBase::Setup ( ) {
  767. ////////////////////////////////////
  768. // First set up grid stuff
  769. if (Cvec) delete [] Cvec;
  770. //if (CSolver) delete [] CSolver;
  771. Cvec = new Eigen::SparseMatrix<Complex, Eigen::RowMajor> [Omegas.size()];
  772. //#ifdef LEMMAUSEOMP
  773. //#pragma omp parallel for schedule(static, 1)
  774. //#endif
  775. for (int iw=0; iw<Omegas.size(); ++iw) {
  776. std::cout << "Build C " << std::endl;
  777. BuildC(sigma, sigma, sigma, iw);
  778. // std::cout << "Build Re(C) " << std::endl;
  779. // BuildCReal(sigma, sigma, sigma, iw);
  780. }
  781. BuildCDirectSolver( ); // templated specialisation
  782. BuildD();
  783. return ;
  784. } // ----- end of method EMSchur3DBase::Setup -----
  785. //--------------------------------------------------------------------------------------
  786. // Class: EMSchur3DBase
  787. // Method: SetLayeredEarthEM
  788. //--------------------------------------------------------------------------------------
  789. void EMSchur3DBase::SetLayeredEarthEM ( std::shared_ptr<LayeredEarthEM> LayModelin ) {
  790. LayModel = LayModelin;
  791. return ;
  792. } // ----- end of method EMSchur3DBase::SetLayeredEarthEM -----
  793. //--------------------------------------------------------------------------------------
  794. // Class: EMSchur3DBase
  795. // Method: PrimaryField
  796. //--------------------------------------------------------------------------------------
  797. void EMSchur3DBase::PrimaryField ( std::shared_ptr<DipoleSource> source, std::shared_ptr<FieldPoints> dpoint ) {
  798. auto EmEarth = Lemma::EMEarth1D::NewSP();
  799. EmEarth->AttachDipoleSource(source);
  800. EmEarth->AttachLayeredEarthEM(LayModel);
  801. EmEarth->AttachFieldPoints(dpoint);
  802. EmEarth->SetFieldsToCalculate(Lemma::E);
  803. EmEarth->SetHankelTransformMethod(ANDERSON801);
  804. std::cout << "Primary Field Calculation" << std::endl;
  805. EmEarth->MakeCalc3();
  806. return ;
  807. } // ----- end of method EMSchur3DBase::PrimaryField -----
  808. //--------------------------------------------------------------------------------------
  809. // Class: EMSchur3DBase
  810. // Method: FillPoints
  811. //--------------------------------------------------------------------------------------
  812. void EMSchur3DBase::FillPoints ( std::shared_ptr<FieldPoints> dpoint) {
  813. dpoint->SetNumberOfPoints(unx + uny + unz);
  814. int ip(0);
  815. Real tx, ty, tz;
  816. Real ox = Grid->GetOx();
  817. Real oy = Grid->GetOy();
  818. Real oz = Grid->GetOz();
  819. tz = oz;
  820. for (int iz=0; iz<nz; ++iz) {
  821. ty = oy;
  822. for (int iy=0; iy<ny; ++iy) {
  823. tx = ox - hx(0)/2.;
  824. for (int ix=0; ix<nx+1; ++ix) {
  825. dpoint->SetLocation(ip, tx, ty, tz);
  826. if (ix < nx) tx += hx(ix);
  827. ++ip;
  828. }
  829. if (iy<ny-1) ty += .5*hy(iy) + .5*hy(iy+1);
  830. }
  831. if (iz<nz-1) tz += .5*hz(iz) + .5*hz(iz+1);
  832. }
  833. tz = oz;
  834. for (int iz=0; iz<nz; ++iz) {
  835. ty = oy - hy(0)/2.;
  836. for (int iy=0; iy<ny+1; ++iy) {
  837. tx = ox;
  838. for (int ix=0; ix<nx; ++ix) {
  839. dpoint->SetLocation(ip, tx, ty, tz);
  840. if (ix<nx-1) tx += .5*hx(ix) + .5*hx(ix+1);
  841. ++ip;
  842. }
  843. if (iy < ny) ty += hy(iy);
  844. }
  845. if (iz<nz-1) tz += .5*hz(iz) + .5*hz(iz+1);
  846. }
  847. tz = oz - hz(0)/2.;
  848. for (int iz=0; iz<nz+1; ++iz) {
  849. ty = oy;
  850. for (int iy=0; iy<ny; ++iy) {
  851. tx = ox;
  852. for (int ix=0; ix<nx; ++ix) {
  853. dpoint->SetLocation(ip, tx, ty, tz);
  854. if (ix<nx-1) tx += .5*hx(ix) + .5*hx(ix+1);
  855. ++ip;
  856. }
  857. if (iy<ny-1) ty += .5*hy(iy) + .5*hy(iy+1);
  858. }
  859. if (iz < nz) tz += hz(iz);
  860. }
  861. return ;
  862. } // ----- end of method EMSchur3DBase::FillPoints -----
  863. //--------------------------------------------------------------------------------------
  864. // Class: EMSchur3DBase
  865. // Method: FillSourceTerms
  866. //--------------------------------------------------------------------------------------
  867. void EMSchur3DBase::FillSourceTerms ( Eigen::Ref<VectorXcr> ms,
  868. Eigen::Ref<VectorXcr> Se, Eigen::Ref<VectorXcr> E0,
  869. std::shared_ptr<FieldPoints> dpoint, const Real& omega ) {
  870. //Complex hsig(0);
  871. //Complex hsigp(0);
  872. Real hsig(0);
  873. Real hsigp(0);
  874. Real EPS(0); //1e-24);
  875. Real EPS2(1e-18);
  876. int iv = 0;
  877. for (int iz=0; iz<nz; ++iz) {
  878. for (int iy=0; iy<ny; ++iy) {
  879. for (int ix=0; ix<nx+1; ++ix) {
  880. // Harmonically average sigma
  881. if (ix == 0) {
  882. hsig = sigma[ix][iy][iz];
  883. hsigp = sigmap[ix][iy][iz];
  884. } else if (ix == nx) {
  885. hsig = sigma[ix-1][iy][iz];
  886. hsigp = sigmap[ix-1][iy][iz];
  887. } else {
  888. hsig = ((hx[ix]+hx[ix-1])/2.) *
  889. (1. / ( (hx[ix ] / (2.*sigma[ix ][iy][iz] + EPS)) +
  890. ( hx[ix-1] / (2.*sigma[ix-1][iy][iz] + EPS)) ) );
  891. hsigp = ((hx[ix]+hx[ix-1])/2.) *
  892. (1. / ( (hx[ix ] / (2.*sigmap[ix ][iy][iz] + EPS)) +
  893. ( hx[ix-1] / (2.*sigmap[ix-1][iy][iz] + EPS)) ) );
  894. }
  895. if (std::abs(hsig) < EPS2) hsig = 0;
  896. if (std::abs(hsigp) < EPS2) hsigp = 0;
  897. //hsig += Complex(0, omega*EPSILON0);
  898. //ioms(iv) = Complex(0, omega*MU0*hsig);
  899. //Se(iv) = ioms(iv)*dpoint->GetEfield(0,iv)(0);
  900. ms(iv) = MU0*hsig;
  901. Se(iv) = MU0*hsigp*
  902. dpoint->GetEfield(0,iv)(0);
  903. E0(iv) = dpoint->GetEfield(0,iv)(0);
  904. ++iv;
  905. }
  906. }
  907. }
  908. for (int iz=0; iz<nz; ++iz) {
  909. for (int iy=0; iy<ny+1; ++iy) {
  910. for (int ix=0; ix<nx; ++ix) {
  911. // Harmonically average sigma
  912. if (iy == 0) {
  913. hsig = sigma[ix][iy][iz];
  914. hsigp = sigmap[ix][iy][iz];
  915. } else if (iy == ny) {
  916. hsig = sigma[ix][iy-1][iz];
  917. hsigp = sigmap[ix][iy-1][iz];
  918. } else {
  919. hsig = ((hy[iy]+hy[iy-1])/2.) *
  920. (1. / ( (hy[iy ] / (2.*sigma[ix][iy ][iz] + EPS)) +
  921. ( hy[iy-1] / (2.*sigma[ix][iy-1][iz] + EPS)) ) );
  922. hsigp = ((hy[iy]+hy[iy-1])/2.) *
  923. (1. / ( (hy[iy ] / (2.*sigmap[ix][iy ][iz] + EPS)) +
  924. ( hy[iy-1] / (2.*sigmap[ix][iy-1][iz] + EPS)) ) );
  925. }
  926. if (std::abs(hsig) < EPS2) hsig = 0;
  927. if (std::abs(hsigp) < EPS2) hsigp = 0;
  928. //hsig += Complex(0, omega*EPSILON0);
  929. //ioms(iv) = Complex(0, omega*MU0*hsig);
  930. //Se(iv) = ioms(iv) * dpoint->GetEfield(0,iv)(1);
  931. ms(iv) = MU0*hsig;
  932. Se(iv) = MU0*hsigp* // TODO was hsigp !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  933. dpoint->GetEfield(0,iv)(1);
  934. E0(iv) = dpoint->GetEfield(0,iv)(1);
  935. ++iv;
  936. }
  937. }
  938. }
  939. for (int iz=0; iz<nz+1; ++iz) {
  940. for (int iy=0; iy<ny; ++iy) {
  941. for (int ix=0; ix<nx; ++ix) {
  942. // Harmonically average sigma
  943. if (iz==0) {
  944. hsig = sigma[ix][iy][iz];
  945. hsigp = sigma[ix][iy][iz];
  946. } else if (iz == nz) {
  947. hsig = sigma[ix][iy][nz-1];
  948. hsigp = sigmap[ix][iy][nz-1];
  949. } else {
  950. hsig = ((hz[iz]+hz[iz-1])/2.) *
  951. (1. / ( (hz[iz ] / (2.*sigma[ix][iy][iz ] + EPS)) +
  952. ( hz[iz-1] / (2.*sigma[ix][iy][iz-1] + EPS)) ) );
  953. hsigp = ((hz[iz]+hz[iz-1])/2.) *
  954. (1. / ( (hz[iz ] / (2.*sigmap[ix][iy][iz ] + EPS)) +
  955. ( hz[iz-1] / (2.*sigmap[ix][iy][iz-1] + EPS)) ) );
  956. }
  957. if (std::abs(hsig) < EPS2) hsig = 0;
  958. if (std::abs(hsigp) < EPS2) hsigp = 0;
  959. //hsig += Complex(0, omega*EPSILON0);
  960. //ioms(iv) = Complex(0, omega*MU0*hsig);
  961. //Se(iv) = ioms(iv) * dpoint->GetEfield(0, iv)(2);
  962. ms(iv) = MU0*hsig;
  963. Se(iv) = MU0*hsig
  964. *dpoint->GetEfield(0, iv)(2);
  965. E0(iv) = dpoint->GetEfield(0, iv)(2);
  966. ++iv;
  967. }
  968. }
  969. }
  970. return ;
  971. } // ----- end of method EMSchur3DBase::FillSourceTerms -----
  972. VectorXcr EMSchur3DBase::StaggeredGridCurl( Eigen::Ref<VectorXcr> A ) {
  973. VectorXcr B = VectorXcr::Zero(3*nx*ny*nz);
  974. VectorXcr dzdy = VectorXcr::Zero( nx*ny*nz);
  975. VectorXcr dydz = VectorXcr::Zero( nx*ny*nz);
  976. VectorXcr dxdz = VectorXcr::Zero( nx*ny*nz);
  977. VectorXcr dzdx = VectorXcr::Zero( nx*ny*nz);
  978. VectorXcr dydx = VectorXcr::Zero( nx*ny*nz);
  979. VectorXcr dxdy = VectorXcr::Zero( nx*ny*nz);
  980. /*
  981. Curl_x = dzdy - dydz
  982. Curl_y = dxdz - dzdx
  983. Curl_z = dydx - dxdy
  984. */
  985. ////////////////////////////////
  986. // X component
  987. int ii = 0;
  988. int iix = 0;
  989. for (int iz=0; iz<nz; ++iz) {
  990. for (int iy=0; iy<ny; ++iy) {
  991. for (int ix=0; ix<nx; ++ix) {
  992. dxdz(iix) = ( A(ii+(nx+1)*ny) - A(ii) ) / hz[iz];
  993. dxdy(iix) = ( A(ii+(nx+1) ) - A(ii) ) / hy[iy];
  994. ++ii;
  995. ++iix;
  996. }
  997. // Jump around boundary
  998. ++ii;
  999. }
  1000. }
  1001. ///////////////////////////////
  1002. // Y component
  1003. int iiy = 0;
  1004. for (int iz=0; iz<nz; ++iz) {
  1005. for (int iy=0; iy<ny; ++iy) {
  1006. for (int ix=0; ix<nx; ++ix) {
  1007. dydz(iiy) = ( A(ii+(nx)*(ny+1)) - A(ii) ) / hz[iz];
  1008. dydx(iiy) = ( A(ii+1 ) - A(ii) ) / hx[ix];
  1009. ++iiy;
  1010. ++ii;
  1011. }
  1012. }
  1013. // Jump around boundary
  1014. ii += nx;
  1015. }
  1016. ////////////////////////////
  1017. // Z component
  1018. int iiz = 0;
  1019. for (int iz=0; iz<nz; ++iz) {
  1020. for (int iy=0; iy<ny; ++iy) {
  1021. for (int ix=0; ix<nx; ++ix) {
  1022. //Az(iiz) = 0.5 * ( A(ii) + A(ii+nx*ny));
  1023. dzdy(iiz) = ( A(ii+nx ) - A(ii) ) / hy[iy];
  1024. dzdx(iiz) = ( A(ii+1 ) - A(ii) ) / hx[ix];
  1025. ++iiz;
  1026. ++ii;
  1027. }
  1028. //++iiz;
  1029. //ii += nx*ny;
  1030. }
  1031. }
  1032. B.segment( 0, nx*ny*nz) = dzdy - dydz;
  1033. B.segment( nx*ny*nz, nx*ny*nz) = dxdz - dzdx;
  1034. B.segment(2*nx*ny*nz, nx*ny*nz) = dydx - dxdy;
  1035. return B;
  1036. }
  1037. //--------------------------------------------------------------------------------------
  1038. // Class: EMSchur3DBase
  1039. // Method: WriteVTKResults
  1040. //--------------------------------------------------------------------------------------
  1041. void EMSchur3DBase::WriteVTKResults ( const std::string& fname, Eigen::Ref<VectorXcr> A,
  1042. Eigen::Ref<VectorXcr> Se, Eigen::Ref<VectorXcr> E0, Eigen::Ref<VectorXcr> E,
  1043. Eigen::Ref<VectorXcr> Phi, Eigen::Ref<VectorXcr> ADiv, Eigen::Ref<VectorXcr> ADiv2,
  1044. Eigen::Ref<VectorXcr> B ) {
  1045. auto VTKGridExporter = RectilinearGridVTKExporter::NewSP();
  1046. VTKGridExporter->SetGrid(Grid);
  1047. // VTK arrays to fill
  1048. // scalars
  1049. vtkDoubleArray *sigmaArray = vtkDoubleArray::New();
  1050. vtkDoubleArray *sigmapArray = vtkDoubleArray::New();
  1051. vtkDoubleArray *phiRealArray = vtkDoubleArray::New();
  1052. vtkDoubleArray *phiImagArray = vtkDoubleArray::New();
  1053. vtkDoubleArray *ADivRealArray = vtkDoubleArray::New();
  1054. vtkDoubleArray *ADivImagArray = vtkDoubleArray::New();
  1055. vtkDoubleArray *ADiv2RealArray = vtkDoubleArray::New();
  1056. vtkDoubleArray *ADiv2ImagArray = vtkDoubleArray::New();
  1057. // vectors
  1058. vtkDoubleArray *AReal = vtkDoubleArray::New();
  1059. vtkDoubleArray *AImag = vtkDoubleArray::New();
  1060. vtkDoubleArray *BReal = vtkDoubleArray::New();
  1061. vtkDoubleArray *BImag = vtkDoubleArray::New();
  1062. vtkDoubleArray *SeReal = vtkDoubleArray::New();
  1063. vtkDoubleArray *SeImag = vtkDoubleArray::New();
  1064. vtkDoubleArray *E0Real = vtkDoubleArray::New();
  1065. vtkDoubleArray *E0Imag = vtkDoubleArray::New();
  1066. vtkDoubleArray *EReal = vtkDoubleArray::New();
  1067. vtkDoubleArray *EImag = vtkDoubleArray::New();
  1068. AReal->SetNumberOfComponents(3);
  1069. AImag->SetNumberOfComponents(3);
  1070. BReal->SetNumberOfComponents(3);
  1071. BImag->SetNumberOfComponents(3);
  1072. SeReal->SetNumberOfComponents(3);
  1073. SeImag->SetNumberOfComponents(3);
  1074. E0Real->SetNumberOfComponents(3);
  1075. E0Imag->SetNumberOfComponents(3);
  1076. EReal->SetNumberOfComponents(3);
  1077. EImag->SetNumberOfComponents(3);
  1078. // Interpolate
  1079. Eigen::VectorXcd Ax(nx*ny*nz); // A
  1080. Eigen::VectorXcd Ay(nx*ny*nz);
  1081. Eigen::VectorXcd Az(nx*ny*nz);
  1082. Eigen::VectorXcd Bx(nx*ny*nz); // B
  1083. Eigen::VectorXcd By(nx*ny*nz);
  1084. Eigen::VectorXcd Bz(nx*ny*nz);
  1085. Eigen::VectorXcd Sex(nx*ny*nz); // Se
  1086. Eigen::VectorXcd Sey(nx*ny*nz);
  1087. Eigen::VectorXcd Sez(nx*ny*nz);
  1088. Eigen::VectorXcd E0x(nx*ny*nz); // E0
  1089. Eigen::VectorXcd E0y(nx*ny*nz);
  1090. Eigen::VectorXcd E0z(nx*ny*nz);
  1091. Eigen::VectorXcd Ex(nx*ny*nz); // E
  1092. Eigen::VectorXcd Ey(nx*ny*nz);
  1093. Eigen::VectorXcd Ez(nx*ny*nz);
  1094. ////////////////////////////////
  1095. // X component
  1096. int ii = 0;
  1097. int iix = 0;
  1098. int ic = 0;
  1099. for (int iz=0; iz<nz; ++iz) {
  1100. for (int iy=0; iy<ny; ++iy) {
  1101. for (int ix=0; ix<nx; ++ix) {
  1102. Ax(iix) = 0.5 * ( A(ii) + A(ii+1));
  1103. Sex(iix) = 0.5 * (Se(ii) + Se(ii+1));
  1104. E0x(iix) = 0.5 * (E0(ii) + E0(ii+1));
  1105. Ex(iix) = 0.5 * ( E(ii) + E(ii+1));
  1106. Bx(iix) = B(ic);
  1107. ++ii;
  1108. ++iix;
  1109. ++ic;
  1110. }
  1111. // Jump around boundary
  1112. ++ii;
  1113. }
  1114. }
  1115. ///////////////////////////////
  1116. // Y component
  1117. int iiy = 0;
  1118. for (int iz=0; iz<nz; ++iz) {
  1119. for (int iy=0; iy<ny; ++iy) {
  1120. for (int ix=0; ix<nx; ++ix) {
  1121. Ay(iiy) = 0.5 * ( A(ii) + A(ii+nx));
  1122. Sey(iiy) = 0.5 * (Se(ii) + Se(ii+nx));
  1123. E0y(iiy) = 0.5 * (E0(ii) + E0(ii+nx));
  1124. Ey(iiy) = 0.5 * ( E(ii) + E(ii+nx));
  1125. By(iiy) = B(ic);
  1126. ++iiy;
  1127. ++ii;
  1128. ++ic;
  1129. }
  1130. }
  1131. // Jump around boundary
  1132. ii += nx;
  1133. }
  1134. ////////////////////////////
  1135. // Z component
  1136. int iiz = 0;
  1137. for (int iz=0; iz<nz; ++iz) {
  1138. for (int iy=0; iy<ny; ++iy) {
  1139. for (int ix=0; ix<nx; ++ix) {
  1140. Az(iiz) = 0.5 * ( A(ii) + A(ii+nx*ny));
  1141. Sez(iiz) = 0.5 * (Se(ii) + Se(ii+nx*ny));
  1142. E0z(iiz) = 0.5 * (E0(ii) + E0(ii+nx*ny));
  1143. Ez(iiz) = 0.5 * ( E(ii) + E(ii+nx*ny));
  1144. Bz(iiz) = B(ic);
  1145. ++iiz;
  1146. ++ii;
  1147. ++ic;
  1148. }
  1149. //++iiz;
  1150. //ii += nx*ny;
  1151. }
  1152. }
  1153. //VTKGridExporter->GetVTKGrid();
  1154. int i = 0;
  1155. for (int iz=0; iz<nz; ++iz) {
  1156. for (int iy=0; iy<ny; ++iy) {
  1157. for (int ix=0; ix<nx; ++ix) {
  1158. // scalars
  1159. //sigmaArray->InsertTuple1(i, std::real(sigma[ix][iy][iz]) );
  1160. //sigmapArray->InsertTuple1(i, std::real(sigmap[ix][iy][iz]) );
  1161. sigmaArray->InsertTuple1(i, sigma[ix][iy][iz] );
  1162. sigmapArray->InsertTuple1(i, sigmap[ix][iy][iz] );
  1163. phiRealArray->InsertTuple1(i, std::real(Phi(i)));
  1164. phiImagArray->InsertTuple1(i, std::imag(Phi(i)));
  1165. ADivRealArray->InsertTuple1(i, std::real(ADiv(i)));
  1166. ADivImagArray->InsertTuple1(i, std::imag(ADiv(i)));
  1167. ADiv2RealArray->InsertTuple1(i, std::real(ADiv2(i)));
  1168. ADiv2ImagArray->InsertTuple1(i, std::imag(ADiv2(i)));
  1169. // vectors
  1170. AReal->InsertTuple3(i, std::real(Ax(i)), std::real(Ay(i)), std::real(Az(i)));
  1171. AImag->InsertTuple3(i, std::imag(Ax(i)), std::imag(Ay(i)), std::imag(Az(i)));
  1172. BReal->InsertTuple3(i, std::real(Bx(i)), std::real(By(i)), std::real(Bz(i)));
  1173. BImag->InsertTuple3(i, std::imag(Bx(i)), std::imag(By(i)), std::imag(Bz(i)));
  1174. SeReal->InsertTuple3(i, std::real(Sex(i)), std::real(Sey(i)), std::real(Sez(i)));
  1175. SeImag->InsertTuple3(i, std::imag(Sex(i)), std::imag(Sey(i)), std::imag(Sez(i)));
  1176. E0Real->InsertTuple3(i, std::real(E0x(i)), std::real(E0y(i)), std::real(E0z(i)));
  1177. E0Imag->InsertTuple3(i, std::imag(E0x(i)), std::imag(E0y(i)), std::imag(E0z(i)));
  1178. EReal-> InsertTuple3(i, std::real(Ex(i)), std::real(Ey(i)), std::real(Ez(i)));
  1179. EImag-> InsertTuple3(i, std::imag(Ex(i)), std::imag(Ey(i)), std::imag(Ez(i)));
  1180. ++i;
  1181. }
  1182. }
  1183. }
  1184. AReal->SetName("A_real");
  1185. AImag->SetName("A_imag");
  1186. BReal->SetName("B_real");
  1187. BImag->SetName("B_imag");
  1188. SeReal->SetName("Se_real");
  1189. SeImag->SetName("Se_imag");
  1190. E0Real->SetName("E0_real");
  1191. E0Imag->SetName("E0_imag");
  1192. EReal->SetName("E_real");
  1193. EImag->SetName("E_imag");
  1194. phiRealArray->SetName("phi Real");
  1195. phiImagArray->SetName("phi imag");
  1196. ADivRealArray->SetName("ini div A real");
  1197. ADivImagArray->SetName("ini div A imag");
  1198. ADiv2RealArray->SetName("div A real");
  1199. ADiv2ImagArray->SetName("div A imag");
  1200. // Add scalar fields
  1201. sigmaArray->SetName("sigma");
  1202. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(sigmaArray);
  1203. sigmapArray->SetName("sigma prime");
  1204. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(sigmapArray);
  1205. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(phiRealArray);
  1206. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(phiImagArray);
  1207. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(ADivRealArray);
  1208. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(ADivImagArray);
  1209. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(ADiv2RealArray);
  1210. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(ADiv2ImagArray);
  1211. // Add Vector Arrays
  1212. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(AReal);
  1213. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(AImag);
  1214. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(SeReal);
  1215. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(SeImag);
  1216. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(E0Real);
  1217. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(E0Imag);
  1218. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(EReal);
  1219. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(EImag);
  1220. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(BReal);
  1221. VTKGridExporter->GetVTKGrid()->GetCellData()->AddArray(BImag);
  1222. // Write
  1223. vtkXMLRectilinearGridWriter *gridWrite = vtkXMLRectilinearGridWriter::New();
  1224. gridWrite->SetInputData( VTKGridExporter->GetVTKGrid() );
  1225. gridWrite->SetFileName( (fname + std::string(".vtr")).c_str());
  1226. gridWrite->Update();
  1227. gridWrite->Write();
  1228. sigmaArray->Delete();
  1229. gridWrite->Delete();
  1230. phiRealArray->Delete();
  1231. phiImagArray->Delete();
  1232. ADivRealArray->Delete();
  1233. ADivImagArray->Delete();
  1234. ADiv2RealArray->Delete();
  1235. ADiv2ImagArray->Delete();
  1236. // TODO delete vectors too!
  1237. AReal->Delete();
  1238. AImag->Delete();
  1239. BReal->Delete();
  1240. BImag->Delete();
  1241. SeReal->Delete();
  1242. SeImag->Delete();
  1243. E0Real->Delete();
  1244. E0Imag->Delete();
  1245. EReal->Delete();
  1246. EImag->Delete();
  1247. return ;
  1248. } // ----- end of method EMSchur3DBase::WriteVTKResults -----
  1249. //--------------------------------------------------------------------------------------
  1250. // Class: EMSchur3DBase
  1251. // Method: LoadMeshToolsConductivityModel
  1252. //--------------------------------------------------------------------------------------
  1253. void EMSchur3DBase::LoadMeshToolsConductivityModel ( const std::string& fname ) {
  1254. if (sigma) Delete3DScalar(sigma);
  1255. //Allocate3DScalar(sigma, Complex(0));
  1256. Allocate3DScalar(sigma, 0.);
  1257. // TODO, there are smarter ways to do this, since the bacground is 1D. But for now, we just have
  1258. // them the same. Eases logic in harmonic averaging. And it's just in one thread, so not so
  1259. // hateful
  1260. if (!LayModel) {
  1261. std::cerr << "In EMSchur3DBase::LoadMeshToolsConductivityModel, you need to set your bacground 1D model\n";
  1262. exit(EXIT_FAILURE);
  1263. }
  1264. if (sigmap) Delete3DScalar(sigmap);
  1265. //Allocate3DScalar(sigmap, Complex(0));
  1266. Allocate3DScalar(sigmap, 0.);
  1267. auto Parser = ASCIIParser::NewSP();
  1268. Parser->Open(fname);
  1269. std::vector<Real> sigmod = Parser->ReadReals( -1 );
  1270. if ( static_cast<int>(sigmod.size()) != nx*ny*nz) {
  1271. std::cerr << "In EMSchur3DBase::LoadMeshToolsConductivityModel model sizes are not in agreement\n";
  1272. exit(EXIT_FAILURE);
  1273. }
  1274. Parser->Close();
  1275. int ic(0);
  1276. for (int ix=0; ix<nx; ++ix) {
  1277. for (int iy=0; iy<ny; ++iy) {
  1278. Real pz = Grid->GetOz();
  1279. for (int iz=0; iz<nz; ++iz) {
  1280. if ( sigmod[ic] != sigmod[ic] ) {
  1281. std::cout << "sigmod problems " << ic << std::endl;
  1282. }
  1283. sigma[ix][iy][iz] = sigmod[ic];
  1284. //sigmap[ix][iy][iz] = sigmod[ic] - LayModel->GetLayerConductivity(LayModel->GetLayerAtThisDepth(pz));
  1285. sigmap[ix][iy][iz] = sigmod[ic] - LayModel->GetLayerConductivity(LayModel->GetLayerAtThisDepth(pz)).real() ;
  1286. ++ic;
  1287. pz += hz[iz];
  1288. }
  1289. }
  1290. }
  1291. // Marker and cell
  1292. MAC = VectorXi::Zero(nx*ny*nz);
  1293. idx.clear();
  1294. ic = 0;
  1295. for (int iz=0; iz<nz; ++iz) {
  1296. for (int iy=0; iy<ny; ++iy) {
  1297. for (int ix=0; ix<nx; ++ix) {
  1298. if (sigma[ix][iy][iz] > 0.) {
  1299. MAC(ic) = 1;
  1300. idx.push_back(ic);
  1301. }
  1302. ++ic;
  1303. }
  1304. }
  1305. }
  1306. return ;
  1307. }
  1308. //--------------------------------------------------------------------------------------
  1309. // Class: EMSchur3DBase
  1310. // Method: SetAEMSurvey
  1311. //--------------------------------------------------------------------------------------
  1312. void EMSchur3DBase::SetAEMSurvey ( std::shared_ptr<AEMSurvey> SurveyIn ) {
  1313. Survey = SurveyIn;
  1314. // determine how many frequencies we are dealing with
  1315. Omegas = 2.*PI*Survey->GetFrequencies();
  1316. return ;
  1317. } // ----- end of method EMSchur3DBase::SetAEMSurvey -----// ----- end of method EMSchur3DBase::LoadMeshToolsConductivityModel -----
  1318. //--------------------------------------------------------------------------------------
  1319. // Class: EMSchur3DBase
  1320. // Method: Solve
  1321. //--------------------------------------------------------------------------------------
  1322. void EMSchur3DBase::Solve ( ) {
  1323. // TODO, these should throw exceptions for a GUI to catch
  1324. if (!LayModel or !Survey or !Grid or !sigma) {
  1325. std::cerr << "You need to set something. EMSChur3D::Solve()";
  1326. exit(EXIT_FAILURE);
  1327. }
  1328. //BuildD(); // strangely necessary, bug, track down.
  1329. Setup();
  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. return ;
  1338. } // ----- end of method EMSchur3DBase::Solve -----
  1339. } // ----- end of Lemma name -----