Surface NMR forward modelling
Nelze vybrat více než 25 témat Téma musí začínat písmenem nebo číslem, může obsahovat pomlčky („-“) a může být dlouhé až 35 znaků.

KernelV0.cpp 20KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486
  1. /* This file is part of Lemma, a geophysical modelling and inversion API.
  2. * More information is available at http://lemmasoftware.org
  3. */
  4. /* This Source Code Form is subject to the terms of the Mozilla Public
  5. * License, v. 2.0. If a copy of the MPL was not distributed with this
  6. * file, You can obtain one at http://mozilla.org/MPL/2.0/.
  7. */
  8. /**
  9. * @file
  10. * @date 11/11/2016 01:47:25 PM
  11. * @version $Id$
  12. * @author Trevor Irons (ti)
  13. * @email tirons@egi.utah.edu
  14. * @copyright Copyright (c) 2016, University of Utah
  15. * @copyright Copyright (c) 2016, Lemma Software, LLC
  16. */
  17. #include "KernelV0.h"
  18. #include "FieldPoints.h"
  19. namespace Lemma {
  20. // ==================== FRIEND METHODS =====================
  21. std::ostream &operator << (std::ostream &stream, const KernelV0 &ob) {
  22. stream << ob.Serialize() << "\n---\n"; // End of doc ---
  23. return stream;
  24. }
  25. // ==================== LIFECYCLE =======================
  26. //--------------------------------------------------------------------------------------
  27. // Class: KernelV0
  28. // Method: KernelV0
  29. // Description: constructor (locked)
  30. //--------------------------------------------------------------------------------------
  31. KernelV0::KernelV0 (const ctor_key&) : LemmaObject( ) {
  32. } // ----- end of method KernelV0::KernelV0 (constructor) -----
  33. //--------------------------------------------------------------------------------------
  34. // Class: KernelV0
  35. // Method: KernelV0
  36. // Description: DeSerializing constructor (locked)
  37. //--------------------------------------------------------------------------------------
  38. KernelV0::KernelV0 (const YAML::Node& node, const ctor_key&) : LemmaObject(node) {
  39. } // ----- end of method KernelV0::KernelV0 (constructor) -----
  40. //--------------------------------------------------------------------------------------
  41. // Class: KernelV0
  42. // Method: NewSP()
  43. // Description: public constructor returing a shared_ptr
  44. //--------------------------------------------------------------------------------------
  45. std::shared_ptr< KernelV0 > KernelV0::NewSP() {
  46. return std::make_shared< KernelV0 >( ctor_key() );
  47. }
  48. //--------------------------------------------------------------------------------------
  49. // Class: KernelV0
  50. // Method: ~KernelV0
  51. // Description: destructor (protected)
  52. //--------------------------------------------------------------------------------------
  53. KernelV0::~KernelV0 () {
  54. } // ----- end of method KernelV0::~KernelV0 (destructor) -----
  55. //--------------------------------------------------------------------------------------
  56. // Class: KernelV0
  57. // Method: Serialize
  58. //--------------------------------------------------------------------------------------
  59. YAML::Node KernelV0::Serialize ( ) const {
  60. YAML::Node node = LemmaObject::Serialize();
  61. node.SetTag( GetName() );
  62. // Coils Transmitters & Receivers
  63. for ( auto txm : TxRx) {
  64. node[txm.first] = txm.second->Serialize();
  65. }
  66. // LayeredEarthEM
  67. node["SigmaModel"] = SigmaModel->Serialize();
  68. return node;
  69. } // ----- end of method KernelV0::Serialize -----
  70. //--------------------------------------------------------------------------------------
  71. // Class: KernelV0
  72. // Method: DeSerialize
  73. //--------------------------------------------------------------------------------------
  74. std::shared_ptr<KernelV0> KernelV0::DeSerialize ( const YAML::Node& node ) {
  75. if (node.Tag() != "KernelV0" ) {
  76. throw DeSerializeTypeMismatch( "KernelV0", node.Tag());
  77. }
  78. return std::make_shared< KernelV0 > ( node, ctor_key() );
  79. } // ----- end of method KernelV0::DeSerialize -----
  80. //--------------------------------------------------------------------------------------
  81. // Class: KernelV0
  82. // Method: DeSerialize
  83. //--------------------------------------------------------------------------------------
  84. void KernelV0::CalculateK0 (const std::vector< std::string>& Tx, const std::vector<std::string >& Rx,
  85. bool vtkOutput ) {
  86. // All EM calculations will share same field points
  87. cpoints = FieldPoints::NewSP();
  88. cpoints->SetNumberOfPoints(8);
  89. for (auto tx : Tx) {
  90. // Set up EMEarth
  91. EMEarths[tx] = EMEarth1D::NewSP();
  92. EMEarths[tx]->AttachWireAntenna(TxRx[tx]);
  93. EMEarths[tx]->AttachLayeredEarthEM(SigmaModel);
  94. EMEarths[tx]->AttachFieldPoints( cpoints );
  95. EMEarths[tx]->SetFieldsToCalculate(H);
  96. // TODO query for method, altough with flat antennae, this is fastest
  97. EMEarths[tx]->SetHankelTransformMethod(ANDERSON801);
  98. EMEarths[tx]->SetTxRxMode(TX);
  99. }
  100. for (auto rx : Rx) {
  101. if (EMEarths.count(rx)) {
  102. EMEarths[rx]->SetTxRxMode(TXRX);
  103. } else {
  104. EMEarths[rx] = EMEarth1D::NewSP();
  105. EMEarths[rx]->AttachWireAntenna(TxRx[rx]);
  106. EMEarths[rx]->AttachLayeredEarthEM(SigmaModel);
  107. EMEarths[rx]->AttachFieldPoints( cpoints );
  108. EMEarths[rx]->SetFieldsToCalculate(H);
  109. // TODO query for method, altough with flat antennae, this is fastest
  110. EMEarths[rx]->SetHankelTransformMethod(ANDERSON801);
  111. EMEarths[rx]->SetTxRxMode(RX);
  112. }
  113. }
  114. IntegrateOnOctreeGrid( 1e-7, vtkOutput );
  115. }
  116. //--------------------------------------------------------------------------------------
  117. // Class: KernelV0
  118. // Method: IntegrateOnOctreeGrid
  119. //--------------------------------------------------------------------------------------
  120. void KernelV0::IntegrateOnOctreeGrid( const Real& tolerance, bool vtkOutput) {
  121. this->tol = tolerance;
  122. //Vector3r Size;
  123. Size << 200,200,20;
  124. //Vector3r Origin;
  125. Origin << 0,0,5.0;
  126. Vector3r cpos; // centre position
  127. //cpos << 100,100,50;
  128. cpos = (Size-Origin).array() / 2.;
  129. int maxlevel;
  130. SUM = 0;
  131. VOLSUM = 0;
  132. nleaves = 0;
  133. if (!vtkOutput) {
  134. EvaluateKids( Size, 0, cpos, 1e6 );
  135. } else {
  136. #ifdef LEMMAUSEVTK
  137. vtkHyperOctree* oct = vtkHyperOctree::New();
  138. oct->SetDimension(3);
  139. oct->SetOrigin( Origin(0), Origin(1), Origin(2) );
  140. oct->SetSize( Size(0), Size(1), Size(2) );
  141. vtkHyperOctreeCursor* curse = oct->NewCellCursor();
  142. curse->ToRoot();
  143. EvaluateKids2( Size, 0, cpos, 1e6, oct, curse );
  144. // Fill in leaf data
  145. vtkDoubleArray* kr = vtkDoubleArray::New();
  146. kr->SetNumberOfComponents(1);
  147. kr->SetName("Re($K_0$)");
  148. kr->SetNumberOfTuples( oct->GetNumberOfLeaves() );
  149. vtkDoubleArray* ki = vtkDoubleArray::New();
  150. ki->SetNumberOfComponents(1);
  151. ki->SetName("Im($K_0$)");
  152. ki->SetNumberOfTuples( oct->GetNumberOfLeaves() );
  153. for (auto leaf : LeafDict) {
  154. kr->InsertTuple1( leaf.first, std::real(leaf.second) );
  155. ki->InsertTuple1( leaf.first, std::imag(leaf.second) );
  156. }
  157. oct->GetLeafData()->AddArray(kr);
  158. oct->GetLeafData()->AddArray(ki);
  159. auto write = vtkXMLHyperOctreeWriter::New();
  160. //write.SetDataModeToAscii()
  161. write->SetInputData(oct);
  162. write->SetFileName("octree.vto");
  163. write->Write();
  164. write->Delete();
  165. kr->Delete();
  166. ki->Delete();
  167. curse->Delete();
  168. oct->Delete();
  169. #else
  170. throw std::runtime_error("IntegrateOnOctreeGrid with vtkOutput requires Lemma with VTK support");
  171. #endif
  172. }
  173. std::cout << "\nVOLSUM=" << VOLSUM << "\tActual=" << Size(0)*Size(1)*Size(2)
  174. << "\tDifference=" << VOLSUM - (Size(0)*Size(1)*Size(2)) << std::endl;
  175. std::cout << "nleaves\t" << nleaves << std::endl;
  176. std::cout << "KSUM\t" << SUM << std::endl;
  177. }
  178. //--------------------------------------------------------------------------------------
  179. // Class: KernelV0
  180. // Method: f
  181. //--------------------------------------------------------------------------------------
  182. Complex KernelV0::f( const Vector3r& r, const Real& volume, const Vector3cr& Ht, const Vector3cr& Hr ) {
  183. return Complex(volume*Ht.dot(Hr));
  184. //return ComputeV0Cell(MU0*Ht, MU0*Hr, volume, 1.0);
  185. }
  186. //--------------------------------------------------------------------------------------
  187. // Class: KernelV0
  188. // Method: ComputeV0Cell
  189. //--------------------------------------------------------------------------------------
  190. Complex KernelV0::ComputeV0Cell(const Vector3cr& Bt,
  191. const Vector3cr& Br, const Real& vol, const Real& phi) {
  192. // Compute the elliptic fields
  193. Vector3r B0hat = {1,0,0};
  194. Vector3r B0 = 53000 * B0hat; // nT
  195. EllipticB EBT = EllipticFieldRep(Bt, B0hat);
  196. EllipticB EBR = EllipticFieldRep(Br, B0hat);
  197. // Compute Mn0
  198. Vector3r Mn0 = ComputeMn0(phi, B0);
  199. Real Mn0Abs = Mn0.norm();
  200. Real Taup = 0.020; // s
  201. Real Ip = 10; // A
  202. // Compute the tipping angle
  203. Real sintheta = std::sin(0.5*GAMMA*Ip*Taup*std::abs(EBT.alpha-EBT.beta));
  204. // Compute phase delay, TODO add transmiiter phase and delay time phase!
  205. Real phase = EBR.zeta+EBT.zeta;
  206. return ComputeV0Cell(EBT, EBR, sintheta, phase, Mn0Abs, vol);
  207. }
  208. //--------------------------------------------------------------------------------------
  209. // Class: KernelV0
  210. // Method: ComputeV0Cell
  211. //--------------------------------------------------------------------------------------
  212. Complex KernelV0::ComputeV0Cell(const EllipticB& EBT, const EllipticB& EBR,
  213. const Real& sintheta, const Real& phase, const Real& Mn0Abs,
  214. const Real& vol) {
  215. Vector3r B0hat = {1,0,0};
  216. Real Larmor = 2500.*2.*PI;
  217. // earth response of receiver adjoint field
  218. Complex ejztr = std::exp(Complex(0, EBR.zeta + EBT.zeta));
  219. Complex PhaseTerm = EBR.bhat.dot(EBT.bhat) +
  220. (B0hat.dot(EBR.bhat.cross(EBT.bhat) ));
  221. return -vol*Complex(0,Larmor)*Mn0Abs*(EBR.alpha+EBR.beta)*ejztr*sintheta*PhaseTerm;
  222. }
  223. Vector3r KernelV0::ComputeMn0(const Real& Porosity, const Vector3r& B0) {
  224. Real Temperature = 283; // in K
  225. Real chi_n = NH2O*((GAMMA*GAMMA*HBAR*HBAR)/(4.*KB*Temperature));
  226. return chi_n*Porosity*B0;
  227. }
  228. //--------------------------------------------------------------------------------------
  229. // Class: KernelV0
  230. // Method: ComputeV0Cell
  231. //--------------------------------------------------------------------------------------
  232. EllipticB KernelV0::EllipticFieldRep (const Vector3cr& B, const Vector3r& B0hat) {
  233. EllipticB ElipB = EllipticB();
  234. Vector3cr Bperp = B.array() - B0hat.dot(B)*B0hat.array();
  235. Real BperpNorm = Bperp.norm();
  236. Complex Bp2 = Bperp.transpose() * Bperp;
  237. VectorXcr iB0 = Complex(0,1)*B0hat.cast<Complex>().array();
  238. ElipB.eizt = std::sqrt(Bp2 / std::abs(Bp2));
  239. ElipB.alpha = INVSQRT2*std::sqrt(BperpNorm*BperpNorm + std::abs(Bp2));
  240. ElipB.beta = sgn(std::real(iB0.dot(Bperp.cross(Bperp.conjugate())))) *
  241. (INVSQRT2)*std::sqrt(std::abs(BperpNorm*BperpNorm-std::abs(Bp2)));
  242. ElipB.bhat = ((Real)1./ElipB.alpha)*(((Real)1./ElipB.eizt)*Bperp.array()).real().array();
  243. ElipB.bhatp = B0hat.cross(ElipB.bhat);
  244. ElipB.zeta = std::real(std::log(ElipB.eizt)/Complex(0,1));
  245. return ElipB;
  246. }
  247. //--------------------------------------------------------------------------------------
  248. // Class: KernelV0
  249. // Method: EvaluateKids
  250. //--------------------------------------------------------------------------------------
  251. bool KernelV0::EvaluateKids( const Vector3r& size, const int& level, const Vector3r& cpos,
  252. const Complex& parentVal ) {
  253. std::cout << "\r" << (int)(1e2*VOLSUM/(Size[0]*Size[1]*Size[2])) << "\t" << nleaves;
  254. std::cout.flush();
  255. // Next level step, interested in one level below
  256. // bitshift requires one extra, faster than, and equivalent to std::pow(2, level+1)
  257. Vector3r step = size.array() / (Real)(1 << (level+1) );
  258. Vector3r step2 = size.array() / (Real)(1 << (level+2) );
  259. Real vol = (step2(0)*step2(1)*step2(2)); // volume of each child
  260. Vector3r pos = cpos - step/2.;
  261. Eigen::Matrix<Real, 8, 3> posadd = (Eigen::Matrix<Real, 8, 3>() <<
  262. 0, 0, 0,
  263. step[0], 0, 0,
  264. 0, step[1], 0,
  265. step[0], step[1], 0,
  266. 0, 0, step[2],
  267. step[0], 0, step[2],
  268. 0, step[1], step[2],
  269. step[0], step[1], step[2] ).finished();
  270. VectorXcr kvals(8); // individual kernel vals
  271. cpoints->ClearFields();
  272. for (int ichild=0; ichild<8; ++ichild) {
  273. Vector3r cp = pos; // Eigen complains about combining these
  274. cp += posadd.row(ichild);
  275. cpoints->SetLocation( ichild, cp );
  276. }
  277. Eigen::Matrix<Complex, 3, 8> Ht = Eigen::Matrix<Complex, 3, 8>::Zero();
  278. Eigen::Matrix<Complex, 3, 8> Hr = Eigen::Matrix<Complex, 3, 8>::Zero();
  279. //Eigen::Matrix< Complex, 8, 3 > Bt;
  280. for ( auto EMCalc : EMEarths ) {
  281. //EMCalc->GetFieldPoints()->ClearFields();
  282. EMCalc.second->CalculateWireAntennaFields();
  283. switch (EMCalc.second->GetTxRxMode()) {
  284. case TX:
  285. Ht += EMCalc.second->GetFieldPoints()->GetHfield(0);
  286. break;
  287. case RX:
  288. Hr += EMCalc.second->GetFieldPoints()->GetHfield(0);
  289. break;
  290. case TXRX:
  291. Ht += EMCalc.second->GetFieldPoints()->GetHfield(0);
  292. Hr += EMCalc.second->GetFieldPoints()->GetHfield(0);
  293. break;
  294. default:
  295. break;
  296. }
  297. }
  298. for (int ichild=0; ichild<8; ++ichild) {
  299. Vector3r cp = pos; // Eigen complains about combining these
  300. cp += posadd.row(ichild);
  301. kvals(ichild) = f(cp, vol, Ht.col(ichild), Hr.col(ichild));
  302. }
  303. Complex ksum = kvals.sum(); // Kernel sum
  304. // Evaluate whether or not furthur splitting is needed
  305. if ( std::abs(ksum - parentVal) > tol || level < 2 ) {
  306. for (int ichild=0; ichild<8; ++ichild) {
  307. Vector3r cp = pos; // Eigen complains about combining these
  308. cp += posadd.row(ichild);
  309. bool isleaf = EvaluateKids( size, level+1, cp, kvals(ichild) );
  310. if (isleaf) { // Include result in final integral
  311. SUM += ksum;
  312. VOLSUM += 8.*vol;
  313. nleaves += 1;
  314. }
  315. }
  316. return false; // not leaf
  317. }
  318. // Save here instead?
  319. return true; // leaf
  320. }
  321. #ifdef LEMMAUSEVTK
  322. //--------------------------------------------------------------------------------------
  323. // Class: KernelV0
  324. // Method: EvaluateKids2 -- same as Evaluate Kids, but include VTK octree generation
  325. //--------------------------------------------------------------------------------------
  326. bool KernelV0::EvaluateKids2( const Vector3r& size, const int& level, const Vector3r& cpos,
  327. const Complex& parentVal, vtkHyperOctree* oct, vtkHyperOctreeCursor* curse) {
  328. std::cout << "\r" << (int)(1e2*VOLSUM/(Size[0]*Size[1]*Size[2])) << "\t" << nleaves;
  329. std::cout.flush();
  330. // Next level step, interested in one level below
  331. // bitshift requires one extra, faster than, and equivalent to std::pow(2, level+1)
  332. Vector3r step = size.array() / (Real)(1 << (level+1) );
  333. Vector3r step2 = size.array() / (Real)(1 << (level+2) );
  334. Real vol = (step2(0)*step2(1)*step2(2)); // volume of each child
  335. Vector3r pos = cpos - step/2.;
  336. Eigen::Matrix<Real, 8, 3> posadd = (Eigen::Matrix<Real, 8, 3>() <<
  337. 0, 0, 0,
  338. step[0], 0, 0,
  339. 0, step[1], 0,
  340. step[0], step[1], 0,
  341. 0, 0, step[2],
  342. step[0], 0, step[2],
  343. 0, step[1], step[2],
  344. step[0], step[1], step[2] ).finished();
  345. VectorXcr kvals(8); // individual kernel vals
  346. cpoints->ClearFields();
  347. for (int ichild=0; ichild<8; ++ichild) {
  348. Vector3r cp = pos; // Eigen complains about combining these
  349. cp += posadd.row(ichild);
  350. cpoints->SetLocation( ichild, cp );
  351. }
  352. Eigen::Matrix<Complex, 3, 8> Ht = Eigen::Matrix<Complex, 3, 8>::Zero();
  353. Eigen::Matrix<Complex, 3, 8> Hr = Eigen::Matrix<Complex, 3, 8>::Zero();
  354. for ( auto EMCalc : EMEarths ) {
  355. //EMCalc->GetFieldPoints()->ClearFields();
  356. EMCalc.second->CalculateWireAntennaFields();
  357. switch (EMCalc.second->GetTxRxMode()) {
  358. case TX:
  359. Ht += EMCalc.second->GetFieldPoints()->GetHfield(0);
  360. break;
  361. case RX:
  362. Hr += EMCalc.second->GetFieldPoints()->GetHfield(0);
  363. break;
  364. case TXRX:
  365. Ht += EMCalc.second->GetFieldPoints()->GetHfield(0);
  366. Hr += EMCalc.second->GetFieldPoints()->GetHfield(0);
  367. break;
  368. default:
  369. break;
  370. }
  371. }
  372. for (int ichild=0; ichild<8; ++ichild) {
  373. Vector3r cp = pos; // Eigen complains about combining these
  374. cp += posadd.row(ichild);
  375. kvals(ichild) = f(cp, vol, Ht.col(ichild), Hr.col(ichild));
  376. }
  377. Complex ksum = kvals.sum(); // Kernel sum
  378. // Evaluate whether or not furthur splitting is needed
  379. if ( std::abs(ksum - parentVal) > tol || level < 2 ) {
  380. oct->SubdivideLeaf(curse);
  381. for (int ichild=0; ichild<8; ++ichild) {
  382. curse->ToChild(ichild);
  383. Vector3r cp = pos; // Eigen complains about combining these
  384. cp += posadd.row(ichild);
  385. /* Test for position via alternative means
  386. Real p[3];
  387. GetPosition(curse, p);
  388. std::cout << cp[0] << "\t" << p[0] << "\t" << cp[1] << "\t" << p[1]
  389. << "\t" << cp[2] << "\t" << p[2] << "\t" << vol<< std::endl;
  390. */
  391. bool isleaf = EvaluateKids2( size, level+1, cp, kvals(ichild), oct, curse );
  392. if (isleaf) { // Include result in final integral
  393. LeafDict[curse->GetLeafId()] = kvals(ichild); // VTK
  394. SUM += ksum;
  395. VOLSUM += 8*vol;
  396. nleaves += 1;
  397. }
  398. curse->ToParent();
  399. }
  400. return false; // not leaf
  401. }
  402. return true; // leaf
  403. }
  404. //--------------------------------------------------------------------------------------
  405. // Class: KernelV0
  406. // Method: GetPosition
  407. //--------------------------------------------------------------------------------------
  408. void KernelV0::GetPosition( vtkHyperOctreeCursor* Cursor, Real* p ) {
  409. Real ratio=1.0/(1<<(Cursor->GetCurrentLevel()));
  410. //step = ((Size).array() / std::pow(2.,Cursor->GetCurrentLevel()));
  411. p[0]=(Cursor->GetIndex(0)+.5)*ratio*this->Size[0]+this->Origin[0] ;//+ .5*step[0];
  412. p[1]=(Cursor->GetIndex(1)+.5)*ratio*this->Size[1]+this->Origin[1] ;//+ .5*step[1];
  413. p[2]=(Cursor->GetIndex(2)+.5)*ratio*this->Size[2]+this->Origin[2] ;//+ .5*step[2];
  414. }
  415. #endif
  416. } // ---- end of namespace Lemma ----
  417. /* vim: set tabstop=4 expandtab */
  418. /* vim: set filetype=cpp */