Lemma is an Electromagnetics API
Вы не можете выбрать более 25 тем Темы должны начинаться с буквы или цифры, могут содержать дефисы(-) и должны содержать не более 35 символов.

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332
  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. /* Original code is port of algorithm published by Key2011
  9. %------------------------------------------------------------------%
  10. % Copyright (c) 2012 by the Society of Exploration Geophysicists. %
  11. % For more information, go to http://software.seg.org/2012/0003 . %
  12. % You must read and accept usage terms at: %
  13. % http://software.seg.org/disclaimer.txt before use. %
  14. %------------------------------------------------------------------%
  15. */
  16. /**
  17. * @file
  18. * @date 02/12/2014 10:28:15 AM
  19. * @version $Id$
  20. * @author Trevor Irons (ti)
  21. * @email Trevor.Irons@xri-geo.com
  22. * @copyright Copyright (c) 2014, XRI Geophysics, LLC
  23. * @copyright Copyright (c) 2014, Trevor Irons
  24. */
  25. #include "QWEKey.h"
  26. //#include <Eigen/Eigenvalues>
  27. namespace Lemma {
  28. // ==================== FRIEND METHODS =====================
  29. std::ostream &operator<<(std::ostream &stream, const QWEKey &ob) {
  30. stream << *(HankelTransform*)(&ob);
  31. return stream;
  32. }
  33. // ==================== LIFECYCLE =======================
  34. //--------------------------------------------------------------------------------------
  35. // Class: QWEKey
  36. // Method: QWEKey
  37. // Description: constructor (protected)
  38. //--------------------------------------------------------------------------------------
  39. //
  40. QWEKey::QWEKey (const std::string& name) : HankelTransform(name), RelTol(1e-12), AbsTol(1e-32), nQuad(61), nDelay(1),
  41. //QWEKey::QWEKey (const std::string& name) : HankelTransform(name), RelTol(1e-38), AbsTol(1e-48), nQuad(39), nDelay(5),
  42. nIntervalsMax(40) {
  43. BesselWeights( J0 ); // TODO experiment with zero weight (J0, J1) options, should be static one time method
  44. } // ----- end of method QWEKey::QWEKey (constructor) -----
  45. //--------------------------------------------------------------------------------------
  46. // Class: QWEKey
  47. // Method: New()
  48. // Description: public constructor
  49. //--------------------------------------------------------------------------------------
  50. QWEKey* QWEKey::New() {
  51. QWEKey* Obj = new QWEKey("QWEKey");
  52. Obj->AttachTo(Obj);
  53. return Obj;
  54. }
  55. //--------------------------------------------------------------------------------------
  56. // Class: QWEKey
  57. // Method: ~QWEKey
  58. // Description: destructor (protected)
  59. //--------------------------------------------------------------------------------------
  60. QWEKey::~QWEKey () {
  61. } // ----- end of method QWEKey::~QWEKey (destructor) -----
  62. //--------------------------------------------------------------------------------------
  63. // Class: QWEKey
  64. // Method: Delete
  65. // Description: public destructor
  66. //--------------------------------------------------------------------------------------
  67. void QWEKey::Delete() {
  68. this->DetachFrom(this);
  69. }
  70. //--------------------------------------------------------------------------------------
  71. // Class: QWEKey
  72. // Method: Release
  73. // Description: destructor (protected)
  74. //--------------------------------------------------------------------------------------
  75. void QWEKey::Release() {
  76. delete this;
  77. }
  78. //--------------------------------------------------------------------------------------
  79. // Class: QWEKey
  80. // Method: Zgauss
  81. //--------------------------------------------------------------------------------------
  82. Complex QWEKey::Zgauss ( const int &ikk, const EMMODE &imode,
  83. const int &itype, const Real &rho,
  84. const Real &wavef, KernelEm1DBase *Kernel ) {
  85. return Textrap(Kernel->GetManagerIndex(), Tn(Kernel->GetManagerIndex())) ;
  86. } // ----- end of method QWEKey::Zgauss -----
  87. //--------------------------------------------------------------------------------------
  88. // Class: QWEKey
  89. // Method: ComputeRelated
  90. //--------------------------------------------------------------------------------------
  91. void QWEKey::ComputeRelated ( const Real& rho, KernelEm1DBase* Kernel ) {
  92. return ;
  93. } // ----- end of method QWEKey::ComputeRelated -----
  94. //--------------------------------------------------------------------------------------
  95. // Class: QWEKey
  96. // Method: ComputeRelated
  97. //--------------------------------------------------------------------------------------
  98. void QWEKey::ComputeRelated ( const Real& rho, std::vector< KernelEm1DBase* > KernelVec ) {
  99. return ;
  100. } // ----- end of method QWEKey::ComputeRelated -----
  101. //--------------------------------------------------------------------------------------
  102. // Class: QWEKey
  103. // Method: ComputeRelated
  104. //--------------------------------------------------------------------------------------
  105. void QWEKey::ComputeRelated ( const Real& rho, KernelEM1DManager* KernelManagerIn ) {
  106. KernelManager = KernelManagerIn; // OK becauase this is internal and we know what we are doing
  107. Lambda = Bx.array()/rho;
  108. Intervals = xInt.array()/rho;
  109. int nrel = (int)(KernelManager->GetSTLVector().size());
  110. Zans = Eigen::Matrix<Complex, Eigen::Dynamic, Eigen::Dynamic>::Zero(1, nrel);
  111. QWE(rho);
  112. return ;
  113. } // ----- end of method QWEKey::ComputeRelated -----
  114. //--------------------------------------------------------------------------------------
  115. // Class: QWEKey
  116. // Method: GaussQuadWeights
  117. //--------------------------------------------------------------------------------------
  118. void QWEKey::GaussQuadWeights(const int& N) {
  119. std::cerr<<"QWEKey needs work to remove Boost, etc." << std::endl;
  120. // Below works with older Eigen, need to find problem
  121. // VectorXr Nv = VectorXr::LinSpaced(N-1, 1, N-1);
  122. // VectorXr beta = 0.5 / (1.-(2.*Nv.array()).pow(-2)).sqrt();
  123. // MatrixXr T = MatrixXr::Zero(N,N);
  124. // //std::cerr << "Eigen ERROR BELOW, QWEKey.cpp QWEKey::GaussQuadWeights, COMMENTED OUT ";
  125. // T.bottomLeftCorner(N-1, N-1) = beta.asDiagonal();
  126. // Eigen::SelfAdjointEigenSolver<MatrixXr> eig( T.selfadjointView< Eigen::Lower >() ); // PROBLEM LINE
  127. // GaussAbscissa = eig.eigenvalues();
  128. // GaussWeights = 2.*eig.eigenvectors().row(0).array().pow(2);
  129. }
  130. //--------------------------------------------------------------------------------------
  131. // Class: QWEKey
  132. // Method: BesselWeights
  133. //--------------------------------------------------------------------------------------
  134. void QWEKey::BesselWeights ( const sZeroType& sType ) {
  135. #ifdef HAVEBOOSTSPECIALFUNCTIONS
  136. GaussQuadWeights(nQuad); // TODO should this be moved out of initializer?
  137. std::vector<Real> bz;
  138. xInt = VectorXr(nIntervalsMax+1);
  139. xInt(0) = 1e-20;
  140. switch (sType) {
  141. case J0:
  142. boost::math::cyl_bessel_j_zero(0.0, 1, nIntervalsMax, std::back_inserter(bz));
  143. xInt.tail(nIntervalsMax) = VectorXr::Map(&bz[0], nIntervalsMax);
  144. break;
  145. case J1:
  146. boost::math::cyl_bessel_j_zero(1.0, 1, nIntervalsMax, std::back_inserter(bz));
  147. xInt.tail(nIntervalsMax) = VectorXr::Map(&bz[0], nIntervalsMax);
  148. break;
  149. case NPI:
  150. xInt << 1e-20, VectorXr::LinSpaced(nIntervalsMax, 1, nIntervalsMax).array() * PI;
  151. break;
  152. }
  153. VectorXr dx = ( xInt.tail(nIntervalsMax) - xInt.head(nIntervalsMax) ).array() / 2.;
  154. // x = GaussAbscissa
  155. // dx in every row GaussWeights+1 rows, cols = n
  156. // dx[0] dx[1] ... dx[N] Gw[0] Gw[0] ... ndX
  157. // dx[0] dx[1] ... dx[N] Gw[1]
  158. MatrixXr Bxm = (dx.transpose().replicate(GaussAbscissa.size(), 1)).eval().array() *
  159. ((GaussAbscissa.replicate(1, dx.size()).array() + 1.));
  160. Bxm.array() += xInt.head(Bxm.cols()).transpose().replicate( Bxm.rows(), 1 ).array();
  161. Bx = VectorXr::Map( &Bxm(0,0), Bxm.size() );
  162. BJ0 = VectorXr(Bx.size());
  163. BJ1 = VectorXr(Bx.size());
  164. int iw = 0;
  165. for (int ii=0; ii<Bx.size(); ++ii) {
  166. BJ0(ii) = boost::math::cyl_bessel_j(0, Bx(ii)) * GaussWeights(iw);
  167. BJ1(ii) = boost::math::cyl_bessel_j(1, Bx(ii)) * GaussWeights(iw);
  168. ++iw;
  169. if (iw == GaussWeights.size()) iw = 0;
  170. }
  171. return ;
  172. # else
  173. std::cerr << "QWEKey requires Boost functionalility that is missing\n";
  174. #endif
  175. } // ----- end of method QWEKey::BesselWeights -----
  176. //--------------------------------------------------------------------------------------
  177. // Class: QWEKey
  178. // Method: QWE
  179. //--------------------------------------------------------------------------------------
  180. void QWEKey::QWE ( const Real& rho ) {
  181. // TODO, is -1 needed here?
  182. int nTerms = nIntervalsMax - nDelay;// - 1;
  183. int nrel = (int)(KernelManager->GetSTLVector().size());
  184. // TODO GREMLINS LIVE IN HERE
  185. MatrixXcr prev = Eigen::Matrix<Complex, Eigen::Dynamic, Eigen::Dynamic>::Zero(1, nrel);
  186. for (int i=0; i<nDelay; ++i) {
  187. getEyKernel(i, 0, rho);
  188. prev += Zans;
  189. }
  190. // Some of these are complex
  191. TS = MatrixXcr::Zero(nrel, nTerms);
  192. Tn = VectorXi::Zero(nrel);
  193. Textrap = MatrixXcr::Zero(nrel, nTerms);
  194. TrelErr = MatrixXr::Zero(nrel, nTerms);
  195. TabsErr = MatrixXr::Zero(nrel, nTerms);
  196. VectorXi Converged = VectorXi::Zero(nrel);
  197. // is nTerms right, 1 array shifting
  198. for (int i=nDelay; i<nTerms; ++i) {
  199. int n = i-nDelay;
  200. getEyKernel(i, 0, rho);
  201. for (int j=0; j<nrel; ++j) {
  202. if (!Converged(j)) { //continue;
  203. Tn(j) = n; // order of the expansion
  204. TS(j,n+1) = TS(j, n) + Zans(0, j); // working array for transformation
  205. /* Compute the Shanks transform using the Epsilon algorithm:
  206. Structured after Weniger (1989, p26) */
  207. /* TI - some kind ob bug here, shanks extrapolation doesn't buy much for TEM at least */
  208. Complex aux2(0);
  209. for (int k=n+1; k>0; --k) {
  210. Complex aux1 = aux2;
  211. aux2 = TS(j,k-1);
  212. Complex ddff = TS(j,k) - aux2;
  213. if (std::abs(ddff) < std::numeric_limits<Real>::min() ) {
  214. TS(j,k-1) = std::numeric_limits<Real>::max() ;
  215. } else {
  216. TS(j,k-1) = aux1 + 1./ddff;
  217. }
  218. }
  219. // The extrapolated result plus the prev integration term:
  220. Textrap(j,n) = TS(j, (n-1)%2)+prev(0, j);
  221. //Textrap(j,n) = TS(j, n%2 + 1)+prev(0, j);
  222. // Step 3: Analyze for convergence:
  223. if (n > 1) {
  224. TabsErr(j,n) = std::abs( Textrap(j, n) - Textrap(j, n-1));
  225. TrelErr(j,n) = TabsErr(j, n) / std::abs(Textrap(j, n)) ;
  226. Converged(j) = TrelErr(j,n) < RelTol + AbsTol/std::abs(Textrap(j,n));
  227. }
  228. }
  229. }
  230. if ( Converged.all() == 1 ) break;
  231. }
  232. // Trim up results
  233. // Clean up the T structure arrays? We can't really do this
  234. // because they are fixed size, maybe see how they are used and
  235. // init to zero. If they are only summed we are OK.
  236. /*
  237. for (int j = 0; j<nrel; ++j) {:nKernels
  238. n = Tn(j);
  239. T(j).extrap = T(j).extrap(1:n);
  240. T(j).relErr = T(j).relErr(1:n);
  241. T(j).absErr = T(j).absErr(1:n);
  242. }
  243. */
  244. return ;
  245. } // ----- end of method QWEKey::QWE -----
  246. //--------------------------------------------------------------------------------------
  247. // Class: QWEKey
  248. // Method: getEyKernel
  249. //--------------------------------------------------------------------------------------
  250. void QWEKey::getEyKernel ( const int& i, const int& idx, const Real& rho ) {
  251. int bidx = i*nQuad;
  252. int nrel = (int)(KernelManager->GetSTLVector().size());
  253. Eigen::Matrix<Complex, Eigen::Dynamic, Eigen::Dynamic > Zwork =
  254. Eigen::Matrix<Complex, Eigen::Dynamic, Eigen::Dynamic>::Zero(nQuad, nrel);
  255. for (int ik=0; ik<nQuad; ++ik) {
  256. KernelManager->ComputeReflectionCoeffs( Lambda(bidx+ik), idx, rho );
  257. for (int ir2=0; ir2<nrel; ++ir2) {
  258. // Zwork* needed due to sign convention (e^-jwt) of FT in filter weights
  259. Zwork(ik, ir2) = std::conj(KernelManager->GetSTLVector()[ir2]->RelBesselArg(Lambda(bidx+ik)));
  260. }
  261. }
  262. Real bma = (Intervals(i+1)-Intervals(i))/2;
  263. for (int ir2=0; ir2<nrel; ++ir2) {
  264. if (KernelManager->GetSTLVector()[ir2]->GetBesselOrder() == 0) {
  265. Zans(0, ir2) = bma * Zwork.col(ir2).dot( BJ0.segment(bidx, nQuad) ); // / rho;
  266. } else {
  267. Zans(0, ir2) = bma * Zwork.col(ir2).dot( BJ1.segment(bidx, nQuad) ); // / rho;
  268. }
  269. }
  270. // fcount += nQuad
  271. return ;
  272. } // ----- end of method QWEKey::getEyKernel -----
  273. void QWEKey::TestPrivate(const int& N) {
  274. //GaussQuadWeights(N);
  275. //std::cout << "abscissa\n" << GaussAbscissa << std::endl;
  276. //std::cout << "weights\n" << GaussWeights << std::endl;
  277. BesselWeights( J1 );
  278. //BesselZeros(0, N);
  279. std::cout << std::scientific;
  280. std::cout << "BJ0" << BJ0 << std::endl;
  281. std::cout << "BJ1" << BJ1 << std::endl;
  282. //std::cout << "Bess Zero\n" << xInt << std::endl;
  283. }
  284. } // ----- end of Lemma name -----