Surface NMR forward modelling
您最多选择25个主题 主题必须以字母或数字开头,可以包含连字符 (-),并且长度不得超过35个字符

KernelV0.cpp 24KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549
  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. // Set up
  87. Larmor = SigmaModel->GetMagneticFieldMagnitude()*GAMMA; // in rad 2246.*2.*PI;
  88. // All EM calculations will share same field points
  89. cpoints = FieldPoints::NewSP();
  90. cpoints->SetNumberOfPoints(8);
  91. for (auto tx : Tx) {
  92. // Set up EMEarth
  93. EMEarths[tx] = EMEarth1D::NewSP();
  94. EMEarths[tx]->AttachWireAntenna(TxRx[tx]);
  95. EMEarths[tx]->AttachLayeredEarthEM(SigmaModel);
  96. EMEarths[tx]->AttachFieldPoints( cpoints );
  97. EMEarths[tx]->SetFieldsToCalculate(H);
  98. // TODO query for method, altough with flat antennae, this is fastest
  99. EMEarths[tx]->SetHankelTransformMethod(ANDERSON801);
  100. EMEarths[tx]->SetTxRxMode(TX);
  101. TxRx[tx]->SetCurrent(1.);
  102. }
  103. for (auto rx : Rx) {
  104. if (EMEarths.count(rx)) {
  105. EMEarths[rx]->SetTxRxMode(TXRX);
  106. } else {
  107. EMEarths[rx] = EMEarth1D::NewSP();
  108. EMEarths[rx]->AttachWireAntenna(TxRx[rx]);
  109. EMEarths[rx]->AttachLayeredEarthEM(SigmaModel);
  110. EMEarths[rx]->AttachFieldPoints( cpoints );
  111. EMEarths[rx]->SetFieldsToCalculate(H);
  112. // TODO query for method, altough with flat antennae, this is fastest
  113. EMEarths[rx]->SetHankelTransformMethod(ANDERSON801);
  114. EMEarths[rx]->SetTxRxMode(RX);
  115. TxRx[rx]->SetCurrent(1.);
  116. }
  117. }
  118. std::cout << "Calculating K0 kernel\n";
  119. Kern = MatrixXcr::Zero( Interfaces.size()-1, PulseI.size() );
  120. for (ilay=0; ilay<Interfaces.size()-1; ++ilay) {
  121. std::cout << "Layer " << ilay << std::endl; //<< " q " << iq << std::endl;
  122. Size(2) = Interfaces(ilay+1) - Interfaces(ilay);
  123. Origin(2) = Interfaces(ilay);
  124. //for (int iq=0; iq< PulseI.size(); ++iq) {
  125. // Ip = PulseI(iq);
  126. //Kern(ilay, iq) =
  127. IntegrateOnOctreeGrid( 0, vtkOutput );
  128. //}
  129. }
  130. std::cout << "\rFinished KERNEL\n";
  131. std::cout << "real\n";
  132. std::cout << Kern.real() << std::endl;
  133. std::cout << "imag\n";
  134. std::cout << Kern.imag() << std::endl;
  135. }
  136. //--------------------------------------------------------------------------------------
  137. // Class: KernelV0
  138. // Method: IntegrateOnOctreeGrid
  139. //--------------------------------------------------------------------------------------
  140. Complex KernelV0::IntegrateOnOctreeGrid( const int& iq, bool vtkOutput) {
  141. Vector3r cpos = Origin + Size/2.;
  142. SUM = 0;
  143. VOLSUM = 0;
  144. nleaves = 0;
  145. if (!vtkOutput) {
  146. EvaluateKids( Size, 0, cpos, VectorXcr::Ones(PulseI.size()) );
  147. } else {
  148. #ifdef LEMMAUSEVTK
  149. vtkHyperOctree* oct = vtkHyperOctree::New();
  150. oct->SetDimension(3);
  151. oct->SetOrigin( Origin(0), Origin(1), Origin(2) );
  152. oct->SetSize( Size(0), Size(1), Size(2) );
  153. vtkHyperOctreeCursor* curse = oct->NewCellCursor();
  154. curse->ToRoot();
  155. EvaluateKids2( Size, 0, cpos, 1e6, oct, curse );
  156. // Fill in leaf data
  157. vtkDoubleArray* kr = vtkDoubleArray::New();
  158. kr->SetNumberOfComponents(1);
  159. kr->SetName("Re($K_0$)");
  160. kr->SetNumberOfTuples( oct->GetNumberOfLeaves() );
  161. vtkDoubleArray* ki = vtkDoubleArray::New();
  162. ki->SetNumberOfComponents(1);
  163. ki->SetName("Im($K_0$)");
  164. ki->SetNumberOfTuples( oct->GetNumberOfLeaves() );
  165. vtkDoubleArray* km = vtkDoubleArray::New();
  166. km->SetNumberOfComponents(1);
  167. km->SetName("mod($K_0$)");
  168. km->SetNumberOfTuples( oct->GetNumberOfLeaves() );
  169. vtkIntArray* kid = vtkIntArray::New();
  170. kid->SetNumberOfComponents(1);
  171. kid->SetName("ID");
  172. kid->SetNumberOfTuples( oct->GetNumberOfLeaves() );
  173. vtkIntArray* kerr = vtkIntArray::New();
  174. kerr->SetNumberOfComponents(1);
  175. kerr->SetName("nleaf");
  176. //Real LeafVol(0);
  177. for (auto leaf : LeafDict) {
  178. kr->InsertTuple1( leaf.first, std::real(leaf.second) );
  179. ki->InsertTuple1( leaf.first, std::imag(leaf.second) );
  180. km->InsertTuple1( leaf.first, std::abs(leaf.second) );
  181. kid->InsertTuple1( leaf.first, leaf.first );
  182. //LeafVol += std::real(leaf.second);
  183. }
  184. //std::cout << "\n\nLeafVol=" << LeafVol << std::endl;
  185. for (auto leaf : LeafDictIdx) {
  186. kerr->InsertTuple1( leaf.first, leaf.second );
  187. }
  188. oct->GetLeafData()->AddArray(kr);
  189. oct->GetLeafData()->AddArray(ki);
  190. oct->GetLeafData()->AddArray(km);
  191. oct->GetLeafData()->AddArray(kid);
  192. oct->GetLeafData()->AddArray(kerr);
  193. auto write = vtkXMLHyperOctreeWriter::New();
  194. //write.SetDataModeToAscii()
  195. write->SetInputData(oct);
  196. std::string fname = std::string("octree-") + to_string(ilay)
  197. + std::string("-") + to_string(iq) + std::string(".vto");
  198. write->SetFileName(fname.c_str());
  199. write->Write();
  200. write->Delete();
  201. //kerr->Delete();
  202. kid->Delete();
  203. kr->Delete();
  204. ki->Delete();
  205. km->Delete();
  206. curse->Delete();
  207. oct->Delete();
  208. #else
  209. throw std::runtime_error("IntegrateOnOctreeGrid with vtkOutput requires Lemma with VTK support");
  210. #endif
  211. }
  212. std::cout << "\nVOLSUM=" << VOLSUM << "\tActual=" << Size(0)*Size(1)*Size(2)
  213. << "\tDifference=" << VOLSUM - (Size(0)*Size(1)*Size(2)) << std::endl;
  214. std::cout << "nleaves\t" << nleaves << std::endl;
  215. std::cout << "KSUM\t" << SUM << std::endl;
  216. return SUM;
  217. }
  218. //--------------------------------------------------------------------------------------
  219. // Class: KernelV0
  220. // Method: f
  221. //--------------------------------------------------------------------------------------
  222. VectorXcr KernelV0::f( const Vector3r& r, const Real& volume, const Vector3cr& Ht, const Vector3cr& Hr ) {
  223. // Compute the elliptic fields
  224. Vector3r B0hat = SigmaModel->GetMagneticFieldUnitVector();
  225. Vector3r B0 = SigmaModel->GetMagneticField();
  226. // Elliptic representation
  227. EllipticB EBT = EllipticFieldRep(MU0*Ht, B0hat);
  228. EllipticB EBR = EllipticFieldRep(MU0*Hr, B0hat);
  229. // Compute Mn0
  230. Vector3r Mn0 = ComputeMn0(1.0, B0);
  231. Real Mn0Abs = Mn0.norm();
  232. // Compute phase delay
  233. // TODO add transmiiter phase and delay time phase!
  234. Real phase = EBR.zeta+EBT.zeta;
  235. // Calcuate vector of all responses
  236. VectorXcr F = VectorXcr::Zero( PulseI.size() );
  237. for (int iq=0; iq<PulseI.size(); ++iq) {
  238. // Compute the tipping angle
  239. Real sintheta = std::sin(0.5*GAMMA*PulseI(iq)*Taup*std::abs(EBT.alpha-EBT.beta));
  240. F(iq) = ComputeV0Cell(EBT, EBR, sintheta, phase, Mn0Abs, volume);
  241. }
  242. return F;
  243. }
  244. // //--------------------------------------------------------------------------------------
  245. // // Class: KernelV0
  246. // // Method: ComputeV0Cell
  247. // //--------------------------------------------------------------------------------------
  248. // Complex KernelV0::ComputeV0Cell(const Vector3cr& Bt,
  249. // const Vector3cr& Br, const Real& vol, const Real& phi) {
  250. //
  251. // // Compute the elliptic fields
  252. // Vector3r B0hat = SigmaModel->GetMagneticFieldUnitVector();
  253. // Vector3r B0 = SigmaModel->GetMagneticField();
  254. //
  255. // // Elliptic representation
  256. // EllipticB EBT = EllipticFieldRep(Bt, B0hat);
  257. // EllipticB EBR = EllipticFieldRep(Br, B0hat);
  258. //
  259. // // Compute Mn0
  260. // Vector3r Mn0 = ComputeMn0(phi, B0);
  261. // Real Mn0Abs = Mn0.norm();
  262. //
  263. // // Compute phase delay, TODO add transmiiter phase and delay time phase!
  264. // Real phase = EBR.zeta+EBT.zeta;
  265. //
  266. // Real sintheta = std::sin(0.5*GAMMA*Ip*Taup*std::abs(EBT.alpha-EBT.beta));
  267. // return 0; ComputeV0Cell(EBT, EBR, sintheta, phase, Mn0Abs, vol);
  268. // }
  269. //--------------------------------------------------------------------------------------
  270. // Class: KernelV0
  271. // Method: ComputeV0Cell
  272. //--------------------------------------------------------------------------------------
  273. Complex KernelV0::ComputeV0Cell(const EllipticB& EBT, const EllipticB& EBR,
  274. const Real& sintheta, const Real& phase, const Real& Mn0Abs,
  275. const Real& vol) {
  276. // earth response of receiver adjoint field
  277. Vector3r B0hat = SigmaModel->GetMagneticFieldUnitVector();
  278. Complex ejztr = std::exp(Complex(0, EBR.zeta + EBT.zeta));
  279. Complex PhaseTerm = EBR.bhat.dot(EBT.bhat) + (B0hat.dot(EBR.bhat.cross(EBT.bhat) ));
  280. return -vol*Complex(0,Larmor)*Mn0Abs*(EBR.alpha+EBR.beta)*ejztr*sintheta*PhaseTerm;
  281. }
  282. //--------------------------------------------------------------------------------------
  283. // Class: KernelV0
  284. // Method: ComputeV0Cell
  285. //--------------------------------------------------------------------------------------
  286. Vector3r KernelV0::ComputeMn0(const Real& Porosity, const Vector3r& B0) {
  287. Real chi_n = NH2O*((GAMMA*GAMMA*HBAR*HBAR)/(4.*KB*Temperature));
  288. return chi_n*Porosity*B0;
  289. }
  290. //--------------------------------------------------------------------------------------
  291. // Class: KernelV0
  292. // Method: ComputeV0Cell
  293. //--------------------------------------------------------------------------------------
  294. EllipticB KernelV0::EllipticFieldRep (const Vector3cr& B, const Vector3r& B0hat) {
  295. EllipticB ElipB = EllipticB();
  296. Vector3cr Bperp = B.array() - B0hat.dot(B)*B0hat.array();
  297. Real BperpNorm = Bperp.norm();
  298. Complex Bp2 = Bperp.transpose() * Bperp;
  299. VectorXcr iB0 = Complex(0,1)*B0hat.cast<Complex>().array();
  300. ElipB.eizt = std::sqrt(Bp2 / std::abs(Bp2));
  301. ElipB.alpha = INVSQRT2*std::sqrt(BperpNorm*BperpNorm + std::abs(Bp2));
  302. ElipB.beta = sgn(std::real(iB0.dot(Bperp.cross(Bperp.conjugate())))) *
  303. (INVSQRT2)*std::sqrt(std::abs(BperpNorm*BperpNorm-std::abs(Bp2)));
  304. ElipB.bhat = ((Real)1./ElipB.alpha)*(((Real)1./ElipB.eizt)*Bperp.array()).real().array();
  305. ElipB.bhatp = B0hat.cross(ElipB.bhat);
  306. ElipB.zeta = std::real(std::log(ElipB.eizt)/Complex(0,1));
  307. return ElipB;
  308. }
  309. //--------------------------------------------------------------------------------------
  310. // Class: KernelV0
  311. // Method: EvaluateKids
  312. //--------------------------------------------------------------------------------------
  313. void KernelV0::EvaluateKids( const Vector3r& size, const int& level, const Vector3r& cpos,
  314. const VectorXcr& parentVal ) {
  315. std::cout << "\r" << (int)(1e2*VOLSUM/(Size[0]*Size[1]*Size[2])) << "\t" << nleaves;
  316. std::cout.flush();
  317. // Next level step, interested in one level below
  318. // bitshift requires one extra, faster than, and equivalent to std::pow(2, level+1)
  319. Vector3r step = size.array() / (Real)(1 << (level+1) );
  320. Real vol = (step(0)*step(1)*step(2)); // volume of each child
  321. Vector3r pos = cpos - step/2.;
  322. Eigen::Matrix<Real, 8, 3> posadd = (Eigen::Matrix<Real, 8, 3>() <<
  323. 0, 0, 0,
  324. step[0], 0, 0,
  325. 0, step[1], 0,
  326. step[0], step[1], 0,
  327. 0, 0, step[2],
  328. step[0], 0, step[2],
  329. 0, step[1], step[2],
  330. step[0], step[1], step[2] ).finished();
  331. //VectorXcr kvals(8); // individual kernel vals
  332. MatrixXcr kvals(8, PulseI.size()); // individual kernel vals
  333. cpoints->ClearFields();
  334. for (int ichild=0; ichild<8; ++ichild) {
  335. Vector3r cp = pos; // Eigen complains about combining these
  336. cp += posadd.row(ichild);
  337. cpoints->SetLocation( ichild, cp );
  338. }
  339. Eigen::Matrix<Complex, 3, 8> Ht = Eigen::Matrix<Complex, 3, 8>::Zero();
  340. Eigen::Matrix<Complex, 3, 8> Hr = Eigen::Matrix<Complex, 3, 8>::Zero();
  341. //Eigen::Matrix< Complex, 8, 3 > Bt;
  342. for ( auto EMCalc : EMEarths ) {
  343. //EMCalc->GetFieldPoints()->ClearFields();
  344. EMCalc.second->CalculateWireAntennaFields();
  345. switch (EMCalc.second->GetTxRxMode()) {
  346. case TX:
  347. Ht += EMCalc.second->GetFieldPoints()->GetHfield(0);
  348. break;
  349. case RX:
  350. Hr += EMCalc.second->GetFieldPoints()->GetHfield(0);
  351. break;
  352. case TXRX:
  353. Ht += EMCalc.second->GetFieldPoints()->GetHfield(0);
  354. Hr += EMCalc.second->GetFieldPoints()->GetHfield(0);
  355. break;
  356. default:
  357. break;
  358. }
  359. }
  360. for (int ichild=0; ichild<8; ++ichild) {
  361. Vector3r cp = pos; // Eigen complains about combining these
  362. cp += posadd.row(ichild);
  363. kvals.row(ichild) = f(cp, vol, Ht.col(ichild), Hr.col(ichild));
  364. }
  365. VectorXcr ksum = kvals.colwise().sum(); // Kernel sum
  366. // Evaluate whether or not furthur splitting is needed
  367. if ( ((ksum - parentVal).array().abs() > tol).any() || level < minLevel && level < maxLevel ) {
  368. // Not a leaf dive further in
  369. for (int ichild=0; ichild<8; ++ichild) {
  370. Vector3r cp = pos; // Eigen complains about combining these
  371. cp += posadd.row(ichild);
  372. EvaluateKids( size, level+1, cp, kvals.row(ichild) );
  373. }
  374. return; // not leaf
  375. }
  376. // implicit else, is a leaf
  377. Kern.row(ilay) += ksum;
  378. VOLSUM += 8.*vol;
  379. nleaves += 1;
  380. return; // is leaf
  381. }
  382. #ifdef LEMMAUSEVTK
  383. //--------------------------------------------------------------------------------------
  384. // Class: KernelV0
  385. // Method: EvaluateKids2 -- same as Evaluate Kids, but include VTK octree generation
  386. //--------------------------------------------------------------------------------------
  387. void KernelV0::EvaluateKids2( const Vector3r& size, const int& level, const Vector3r& cpos,
  388. const Complex& parentVal, vtkHyperOctree* oct, vtkHyperOctreeCursor* curse) {
  389. std::cout << "\r" << (int)(1e2*VOLSUM/(Size[0]*Size[1]*Size[2])) << "\t" << nleaves;
  390. std::cout.flush();
  391. // Next level step, interested in one level below
  392. // bitshift requires one extra, faster than, and equivalent to std::pow(2, level+1)
  393. Vector3r step = size.array() / (Real)(1 << (level+1) );
  394. Real vol = (step(0)*step(1)*step(2)); // volume of each child
  395. Vector3r pos = cpos - step/2.;
  396. Eigen::Matrix<Real, 8, 3> posadd = (Eigen::Matrix<Real, 8, 3>() <<
  397. 0, 0, 0,
  398. step[0], 0, 0,
  399. 0, step[1], 0,
  400. step[0], step[1], 0,
  401. 0, 0, step[2],
  402. step[0], 0, step[2],
  403. 0, step[1], step[2],
  404. step[0], step[1], step[2] ).finished();
  405. VectorXcr kvals(8); // individual kernel vals
  406. cpoints->ClearFields();
  407. for (int ichild=0; ichild<8; ++ichild) {
  408. Vector3r cp = pos; // Eigen complains about combining these
  409. cp += posadd.row(ichild);
  410. cpoints->SetLocation( ichild, cp );
  411. }
  412. Eigen::Matrix<Complex, 3, 8> Ht = Eigen::Matrix<Complex, 3, 8>::Zero();
  413. Eigen::Matrix<Complex, 3, 8> Hr = Eigen::Matrix<Complex, 3, 8>::Zero();
  414. for ( auto EMCalc : EMEarths ) {
  415. //EMCalc->GetFieldPoints()->ClearFields();
  416. EMCalc.second->CalculateWireAntennaFields();
  417. switch (EMCalc.second->GetTxRxMode()) {
  418. case TX:
  419. Ht += EMCalc.second->GetFieldPoints()->GetHfield(0);
  420. break;
  421. case RX:
  422. Hr += EMCalc.second->GetFieldPoints()->GetHfield(0);
  423. break;
  424. case TXRX:
  425. Ht += EMCalc.second->GetFieldPoints()->GetHfield(0);
  426. Hr += EMCalc.second->GetFieldPoints()->GetHfield(0);
  427. break;
  428. default:
  429. break;
  430. }
  431. }
  432. for (int ichild=0; ichild<8; ++ichild) {
  433. Vector3r cp = pos; // Eigen complains about combining these
  434. cp += posadd.row(ichild);
  435. kvals(ichild) = f(cp, vol, Ht.col(ichild), Hr.col(ichild))(0);
  436. }
  437. Complex ksum = kvals.sum(); // Kernel sum
  438. // Evaluate whether or not furthur splitting is needed
  439. if ( std::abs(ksum - parentVal) > tol || level < minLevel && level < maxLevel ) {
  440. //if ( std::abs(ksum.real()-parentVal.real())>tol || std::abs(ksum.imag()-parentVal.imag()) >tol ) {
  441. oct->SubdivideLeaf(curse);
  442. for (int ichild=0; ichild<8; ++ichild) {
  443. curse->ToChild(ichild);
  444. Vector3r cp = pos; // Eigen complains about combining these
  445. cp += posadd.row(ichild);
  446. /* Test for position via alternative means */
  447. /*
  448. Real p[3];
  449. GetPosition(curse, p);
  450. if ( (Vector3r(p) - cp).norm() > 1e-8 ) {
  451. std::cout << "ERROR @ nleaves" << nleaves << "\n" << cp[0] << "\t" << p[0] << "\t" << cp[1] << "\t" << p[1]
  452. << "\t" << cp[2] << "\t" << p[2] << "\t" << vol<< std::endl;
  453. throw std::runtime_error("doom");
  454. }
  455. */
  456. /* End of position test */
  457. EvaluateKids2( size, level+1, cp, kvals(ichild), oct, curse );
  458. curse->ToParent();
  459. }
  460. return; // not a leaf
  461. }
  462. LeafDict[curse->GetLeafId()] = ksum/(8.*vol); //kvals(ichild) / vol; // VTK
  463. LeafDictIdx[curse->GetLeafId()] = nleaves; // VTK
  464. SUM += ksum;
  465. VOLSUM += 8*vol;
  466. nleaves += 1;
  467. return; // is a leaf
  468. }
  469. //--------------------------------------------------------------------------------------
  470. // Class: KernelV0
  471. // Method: GetPosition
  472. //--------------------------------------------------------------------------------------
  473. void KernelV0::GetPosition( vtkHyperOctreeCursor* Cursor, Real* p ) {
  474. Real ratio=1.0/(1<<(Cursor->GetCurrentLevel()));
  475. //step = ((Size).array() / std::pow(2.,Cursor->GetCurrentLevel()));
  476. p[0]=(Cursor->GetIndex(0)+.5)*ratio*this->Size[0]+this->Origin[0] ;//+ .5*step[0];
  477. p[1]=(Cursor->GetIndex(1)+.5)*ratio*this->Size[1]+this->Origin[1] ;//+ .5*step[1];
  478. p[2]=(Cursor->GetIndex(2)+.5)*ratio*this->Size[2]+this->Origin[2] ;//+ .5*step[2];
  479. }
  480. #endif
  481. } // ---- end of namespace Lemma ----
  482. /* vim: set tabstop=4 expandtab */
  483. /* vim: set filetype=cpp */