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.

bicgstab.h 26KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936
  1. // ===========================================================================
  2. //
  3. // Filename: bicgstab.h
  4. //
  5. // Description:
  6. //
  7. // Version: 0.0
  8. // Created: 10/27/2009 03:15:06 PM
  9. // Revision: none
  10. // Compiler: g++ (c++)
  11. //
  12. // Author: Trevor Irons (ti)
  13. //
  14. // Organisation: Colorado School of Mines (CSM)
  15. // United States Geological Survey (USGS)
  16. //
  17. // Email: tirons@mines.edu, tirons@usgs.gov
  18. //
  19. // This program is free software: you can redistribute it and/or modify
  20. // it under the terms of the GNU General Public License as published by
  21. // the Free Software Foundation, either version 3 of the License, or
  22. // (at your option) any later version.
  23. //
  24. // This program is distributed in the hope that it will be useful,
  25. // but WITHOUT ANY WARRANTY; without even the implied warranty of
  26. // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
  27. // GNU General Public License for more details.
  28. //
  29. // You should have received a copy of the GNU General Public License
  30. // along with this program. If not, see <http://www.gnu.org/licenses/>.
  31. //
  32. // ===========================================================================
  33. #include <Eigen/Core>
  34. #include <Eigen/Sparse>
  35. #ifdef CHOLMODPRECONDITION
  36. #include <Eigen/CholmodSupport>
  37. #endif // CHOLMODPRECONDITION
  38. //#include "unsupported/Eigen/IterativeSolvers"
  39. //#include <unsupported/Eigen/SuperLUSupport>
  40. #include <iostream>
  41. #include <string>
  42. #include <complex>
  43. #include <fstream>
  44. #include "lemma.h"
  45. #include "timer.h"
  46. using namespace Eigen;
  47. using namespace Lemma;
  48. //typedef Eigen::VectorXcd VectorXcr;
  49. typedef Eigen::SparseMatrix<Complex> SparseMat;
  50. // On Input
  51. // A = Matrix
  52. // B = Right hand side
  53. // X = initial guess, and solution
  54. // maxit = maximum Number of iterations
  55. // tol = error tolerance
  56. // On Output
  57. // X real solution vector
  58. // errorn = Real error norm
  59. int bicgstab(const SparseMat &A, const SparseMat &M, const VectorXcr &b, VectorXcr &x,
  60. int &max_it, Real &tol, Real &errorn, int &iter_done,
  61. const bool& banner = true) {
  62. Complex omega, rho, rho_1, alpha, beta;
  63. Real bnrm2, eps, errmin;
  64. int n, iter; //, istat;
  65. // Determine size of system and init vectors
  66. n = x.size();
  67. VectorXcr r(n);
  68. VectorXcr r_tld(n);
  69. VectorXcr p(n);
  70. VectorXcr v(n);
  71. VectorXcr p_hat(n);
  72. VectorXcr s(n);
  73. VectorXcr s_hat(n);
  74. VectorXcr t(n);
  75. VectorXcr xmin(n);
  76. if (banner) {
  77. std::cout << "Start BiCGStab, memory needed: "
  78. << (sizeof(Complex)*(9+2)*n/(1024.*1024*1024)) << " [Gb]\n";
  79. }
  80. // Initialise
  81. iter_done = 0;
  82. v.setConstant(0.); // not necessary I don't think
  83. t.setConstant(0.);
  84. eps = 1e-100;
  85. bnrm2 = b.norm();
  86. if (bnrm2 == 0) {
  87. x.setConstant(0.0);
  88. errorn = 0;
  89. std::cerr << "Trivial case of Ax = b, where b is 0\n";
  90. return (0);
  91. }
  92. // If there is an initial guess
  93. if ( x.norm() ) {
  94. r = b - A.selfadjointView<Eigen::Upper>()*x;
  95. //r = b - A*x;
  96. } else {
  97. r = b;
  98. }
  99. errorn = r.norm() / bnrm2;
  100. omega = 1.;
  101. r_tld = r;
  102. errmin = 1e30;
  103. // Get down to business
  104. for (iter=0; iter<max_it; ++iter) {
  105. rho = r_tld.dot(r);
  106. if ( abs(rho) < eps) return (0);
  107. if (iter > 0) {
  108. beta = (rho/rho_1) * (alpha/omega);
  109. p = r.array() + beta*(p.array()-omega*v.array()).array();
  110. } else {
  111. p = r;
  112. }
  113. // Use pseudo inverse to get approximate answer
  114. //#pragma omp sections
  115. p_hat = M*p;
  116. //v = A*p_hat; // TODO double check
  117. v = A.selfadjointView<Eigen::Upper>()*p_hat; // TODO double check
  118. alpha = rho / r_tld.dot(v);
  119. s = r.array() - alpha*v.array();
  120. errorn = s.norm()/bnrm2;
  121. if (errorn < tol && iter > 1) {
  122. x.array() += alpha*p_hat.array();
  123. return (0);
  124. }
  125. s_hat = M*s;
  126. t = A.selfadjointView<Eigen::Upper>()*s_hat;
  127. //t = A*s_hat;
  128. omega = t.dot(s) / t.dot(t);
  129. x.array() += alpha*p_hat.array() + omega*s_hat.array();
  130. r = s.array() - omega*t.array();
  131. errorn = r.norm() / bnrm2;
  132. iter_done = iter;
  133. if (errorn < errmin) {
  134. // remember the model with the smallest norm
  135. errmin = errorn;
  136. xmin = x;
  137. }
  138. if ( errorn <= tol ) return (0);
  139. if ( abs(omega) < eps ) return (0);
  140. rho_1 = rho;
  141. }
  142. return (0);
  143. }
  144. template <typename Preconditioner>
  145. bool preconditionedBiCGStab(const SparseMat &A, const Preconditioner &M,
  146. const Ref< VectorXcr const > b,
  147. Ref <VectorXcr > x,
  148. const int &max_it, const Real &tol,
  149. Real &errorn, int &iter_done) {
  150. Complex omega, rho, rho_1, alpha, beta;
  151. Real bnrm2, eps;
  152. int n, iter;
  153. Real tol2 = tol*tol;
  154. // Determine size of system and init vectors
  155. n = x.size();
  156. VectorXcd r(n);
  157. VectorXcd r_tld(n);
  158. VectorXcd p(n);
  159. VectorXcd s(n);
  160. VectorXcd s_hat(n);
  161. VectorXcd p_hat(n);
  162. VectorXcd v = VectorXcr::Zero(n);
  163. VectorXcd t = VectorXcr::Zero(n);
  164. //std::cout << "Start BiCGStab, memory needed: "
  165. // << (sizeof(Complex)*(8+2)*n/(1024.*1024)) / (1024.) << " [Gb]\n";
  166. // Initialise
  167. iter_done = 0;
  168. eps = 1e-100;
  169. bnrm2 = b.squaredNorm();
  170. if (bnrm2 == 0) {
  171. x.setConstant(0.0);
  172. errorn = 0;
  173. std::cerr << "Trivial case of Ax = b, where b is 0\n";
  174. return (false);
  175. }
  176. // If there is an initial guess
  177. if ( x.squaredNorm() ) {
  178. r = b - A.selfadjointView<Eigen::Upper>()*x;
  179. } else {
  180. r = b;
  181. }
  182. errorn = r.squaredNorm() / bnrm2;
  183. omega = 1.;
  184. r_tld = r;
  185. // Get down to business
  186. for (iter=0; iter<max_it; ++iter) {
  187. rho = r_tld.dot(r);
  188. if (abs(rho) < eps) {
  189. std::cerr << "arbitrary orthogonality issue in bicgstab\n";
  190. std::cerr << "consider eigen restarting\n";
  191. return (false);
  192. }
  193. if (iter > 0) {
  194. beta = (rho/rho_1) * (alpha/omega);
  195. p = r + beta*(p-omega*v);
  196. } else {
  197. p = r;
  198. }
  199. p_hat = M.solve(p);
  200. v.noalias() = A.selfadjointView<Eigen::Upper>()*p_hat;
  201. alpha = rho / r_tld.dot(v);
  202. s = r - alpha*v;
  203. errorn = s.squaredNorm()/bnrm2;
  204. if (errorn < tol2 && iter > 1) {
  205. x = x + alpha*p_hat;
  206. errorn = std::sqrt(errorn);
  207. return (true);
  208. }
  209. s_hat = M.solve(s);
  210. t.noalias() = A.selfadjointView<Eigen::Upper>()*s_hat;
  211. omega = t.dot(s) / t.dot(t);
  212. x += alpha*p_hat + omega*s_hat;
  213. r = s - omega*t;
  214. errorn = r.squaredNorm() / bnrm2;
  215. iter_done = iter;
  216. if ( errorn <= tol2 || abs(omega) < eps) {
  217. errorn = std::sqrt(errorn);
  218. return (true);
  219. }
  220. rho_1 = rho;
  221. }
  222. return (false);
  223. }
  224. template <typename Preconditioner>
  225. bool preconditionedSCBiCG(const SparseMat &A, const Preconditioner &M,
  226. const Ref< VectorXcr const > b,
  227. Ref <VectorXcr > x,
  228. const int &max_iter, const Real &tol,
  229. Real &errorn, int &iter_done) {
  230. Real resid;
  231. VectorXcr p, z, q;
  232. Complex alpha, beta, rho, rho_1;
  233. Real normb = b.norm( );
  234. VectorXcr r = b - A*x;
  235. if (normb == 0.0) normb = 1;
  236. if ((resid = r.norm( ) / normb) <= tol) {
  237. errorn = resid;
  238. iter_done = 0;
  239. return 0;
  240. }
  241. for (int i = 1; i <= max_iter; i++) {
  242. z = M.solve(r);
  243. rho = r.dot(z);
  244. if (i == 1) p = z;
  245. else {
  246. beta = rho / rho_1;
  247. p = z + beta * p;
  248. }
  249. q = A*p;
  250. alpha = rho / p.dot(q);
  251. x += alpha * p;
  252. r -= alpha * q;
  253. std::cout << "resid\t" << resid << std::endl;
  254. if ((resid = r.norm( ) / normb) <= tol) {
  255. errorn = resid;
  256. iter_done = i;
  257. return 0;
  258. }
  259. rho_1 = rho;
  260. }
  261. errorn = resid;
  262. return (false);
  263. }
  264. /** \internal Low-level conjugate gradient algorithm
  265. * \param mat The matrix A
  266. * \param rhs The right hand side vector b
  267. * \param x On input and initial solution, on output the computed solution.
  268. * \param precond A preconditioner being able to efficiently solve for an
  269. * approximation of Ax=b (regardless of b)
  270. * \param iters On input the max number of iteration, on output the number of performed iterations.
  271. * \param tol_error On input the tolerance error, on output an estimation of the relative error.
  272. */
  273. template<typename Rhs, typename Dest, typename Preconditioner>
  274. EIGEN_DONT_INLINE
  275. void conjugateGradient(const SparseMat& mat, const Rhs& rhs, Dest& x,
  276. const Preconditioner& precond, int& iters,
  277. typename Dest::RealScalar& tol_error)
  278. {
  279. using std::sqrt;
  280. using std::abs;
  281. typedef typename Dest::RealScalar RealScalar;
  282. typedef typename Dest::Scalar Scalar;
  283. typedef Matrix<Scalar,Dynamic,1> VectorType;
  284. RealScalar tol = tol_error;
  285. int maxIters = iters;
  286. int n = mat.cols();
  287. VectorType residual = rhs - mat.selfadjointView<Eigen::Upper>() * x; //initial residual
  288. RealScalar rhsNorm2 = rhs.squaredNorm();
  289. if(rhsNorm2 == 0)
  290. {
  291. x.setZero();
  292. iters = 0;
  293. tol_error = 0;
  294. return;
  295. }
  296. RealScalar threshold = tol*tol*rhsNorm2;
  297. RealScalar residualNorm2 = residual.squaredNorm();
  298. if (residualNorm2 < threshold)
  299. {
  300. iters = 0;
  301. tol_error = sqrt(residualNorm2 / rhsNorm2);
  302. return;
  303. }
  304. VectorType p(n);
  305. p = precond.solve(residual); //initial search direction
  306. VectorType z(n), tmp(n);
  307. RealScalar absNew = numext::real(residual.dot(p)); // the square of the absolute value of r scaled by invM
  308. int i = 0;
  309. while(i < maxIters)
  310. {
  311. tmp.noalias() = mat.selfadjointView<Eigen::Upper>() * p; // the bottleneck of the algorithm
  312. Scalar alpha = absNew / p.dot(tmp); // the amount we travel on dir
  313. x += alpha * p; // update solution
  314. residual -= alpha * tmp; // update residue
  315. residualNorm2 = residual.squaredNorm();
  316. if(residualNorm2 < threshold)
  317. break;
  318. z = precond.solve(residual); // approximately solve for "A z = residual"
  319. RealScalar absOld = absNew;
  320. absNew = numext::real(residual.dot(z)); // update the absolute value of r
  321. RealScalar beta = absNew / absOld; // calculate the Gram-Schmidt value used to create the new search direction
  322. p = z + beta * p; // update search direction
  323. i++;
  324. }
  325. tol_error = sqrt(residualNorm2 / rhsNorm2);
  326. iters = i;
  327. }
  328. // // Computes implicit
  329. // VectorXcr implicitDCInvBPhi (const SparseMat& D, const SparseMat& C,
  330. // const SparseMat& B, const SparseMat& MC,
  331. // const VectorXcr& Phi, Real& tol,
  332. // int& max_it) {
  333. // int iter_done(0);
  334. // Real errorn(0);
  335. // VectorXcr b = B*Phi;
  336. // VectorXcr y = VectorXcr::Zero(C.rows()) ; // = C^1*b;
  337. // bicgstab(C, MC, b, y, max_it, tol, errorn, iter_done, false);
  338. // //std::cout << "Temp " << errorn << std::endl;
  339. // return D*y;
  340. // }
  341. // Computes implicit
  342. VectorXcr implicitDCInvBPhi (const SparseMat& D, const SparseMat& C,
  343. const VectorXcr& ioms, const SparseMat& MC,
  344. const VectorXcr& Phi, Real& tol,
  345. int& max_it) {
  346. int iter_done(0);
  347. Real errorn(0);
  348. VectorXcr b = (ioms).asDiagonal() * (D.transpose()*Phi);
  349. VectorXcr y = VectorXcr::Zero(C.rows()) ; // = C^1*b;
  350. bicgstab(C, MC, b, y, max_it, tol, errorn, iter_done, false);
  351. //std::cout << "Temp " << errorn << std::endl;
  352. max_it = iter_done;
  353. return D*y;
  354. }
  355. // Computes implicit
  356. template <typename Preconditioner>
  357. VectorXcr implicitDCInvBPhi2 (const SparseMat& D, const SparseMat& C,
  358. const Ref<VectorXcr const> ioms, const Preconditioner& solver,
  359. const Ref<VectorXcr const> Phi, Real& tol,
  360. int& max_it) {
  361. VectorXcr b = (ioms).asDiagonal() * (D.transpose()*Phi);
  362. VectorXcr y = VectorXcr::Zero(C.rows()) ; // = C^1*b;
  363. // Home Made
  364. //int iter_done(0);
  365. //Real errorn(0);
  366. //preconditionedBiCGStab(C, solver, b, y, max_it, tol, errorn, iter_done); //, false); // Jacobi M
  367. //max_it = iter_done;
  368. // Eigen BiCGStab
  369. Eigen::BiCGSTAB<SparseMatrix<Complex> > BiCG;
  370. BiCG.compute( C ); // TODO move this out of this loop!
  371. y = BiCG.solve(b);
  372. max_it = BiCG.iterations();
  373. tol = BiCG.error();
  374. // Direct
  375. /*
  376. std::cout << "Computing LLT" << std::endl;
  377. Eigen::SimplicialLLT<SparseMatrix<Complex>, Eigen::Upper, Eigen::AMDOrdering<int> > LLT;
  378. LLT.compute(C);
  379. max_it = 1;
  380. std::cout << "Computed LLT" << std::endl;
  381. y = LLT.solve(b);
  382. */
  383. return D*y;
  384. }
  385. // Computes implicit
  386. //template <typename Solver>
  387. template < typename Solver >
  388. inline VectorXcr implicitDCInvBPhi3 (const SparseMat& D, const Solver& solver,
  389. const Ref<VectorXcr const> ioms,
  390. const Ref<VectorXcr const> Phi, Real& tol,
  391. int& max_it) {
  392. VectorXcr b = (ioms).asDiagonal() * (D.transpose()*Phi);
  393. VectorXcr y = solver.solve(b);
  394. max_it = 0;
  395. //max_it = solver.iterations();
  396. //errorn = solver.error();
  397. return D*y;
  398. }
  399. // // Simple extraction of indices in idx into reduceed array x1
  400. // void vmap( const Ref<VectorXcr const> x0, Ref<VectorXcr> x1, const std::vector<int>& idx ) {
  401. // for (unsigned int ii=0; ii<idx.size(); ++ii) {
  402. // x1(ii) = x0(idx[ii]);
  403. // }
  404. // }
  405. // Simple extraction of indices in idx into reduceed array x1
  406. VectorXcr vmap( const Ref<VectorXcr const> x0, const std::vector<int>& idx ) {
  407. VectorXcr x1 = VectorXcr::Zero( idx.size() );
  408. for (unsigned int ii=0; ii<idx.size(); ++ii) {
  409. x1(ii) = x0(idx[ii]);
  410. }
  411. return x1;
  412. }
  413. // reverse of above
  414. void ivmap( Ref<VectorXcr > x0, const Ref<VectorXcr const> x1, const std::vector<int>& idx ) {
  415. for (unsigned int ii=0; ii<idx.size(); ++ii) {
  416. x0(idx[ii]) = x1(ii);
  417. }
  418. }
  419. // On Input
  420. // A = Matrix
  421. // B = Right hand side
  422. // X = initial guess, and solution
  423. // maxit = maximum Number of iterations
  424. // tol = error tolerance
  425. // On Output
  426. // X real solution vector
  427. // errorn = Real error norm
  428. template < typename CSolver >
  429. int implicitbicgstab(//const SparseMat& D,
  430. //const SparseMat& C,
  431. const Ref< Eigen::SparseMatrix<Complex> const > D,
  432. const std::vector<int>& idx,
  433. const Ref< VectorXcr const > ioms,
  434. const Ref< VectorXcr const > rhs,
  435. Ref <VectorXcr> phi,
  436. CSolver& solver,
  437. int &max_it, const Real &tol, Real &errorn, int &iter_done, ofstream& logio) {
  438. logio << "using the preconditioned implicit solver" << std::endl;
  439. Complex omega, rho, rho_1, alpha, beta;
  440. Real tol2;
  441. int iter, max_it2, max_it1;
  442. // Look at reduced problem
  443. VectorXcr rhs2 = vmap(rhs, idx);
  444. VectorXcr phi2 = vmap(phi, idx);
  445. // Determine size of system and init vectors
  446. int n = idx.size(); // was phi.size();
  447. VectorXcr r(n);
  448. VectorXcr r_tld(n);
  449. VectorXcr p(n);
  450. VectorXcr s(n);
  451. VectorXcr v = VectorXcr::Zero(n);
  452. VectorXcr t = VectorXcr::Zero(n);
  453. // TODO, refigure for implicit large system
  454. // std::cout << "Start BiCGStab, memory needed: "
  455. // << (sizeof(Complex)*(9+2)*n/(1024.*1024*1024)) << " [Gb]\n";
  456. // Initialise
  457. iter_done = 0;
  458. Real eps = 1e-100;
  459. Real bnrm2 = rhs.norm();
  460. if (bnrm2 == 0) {
  461. phi.setConstant(0.0);
  462. errorn = 0;
  463. std::cerr << "Trivial case of Ax = b, where b is 0\n";
  464. return (0);
  465. }
  466. // If there is an initial guess
  467. if ( phi.norm() ) {
  468. tol2 = tol;
  469. max_it2 = 50000;
  470. //r = rhs - implicitDCInvBPhi3(D, solver, ioms, phi, tol2, max_it2);
  471. //r = rhs - implicitDCInvBPhi3(D, solver, ioms, phi, tol2, max_it2);
  472. r = rhs2 - vmap( implicitDCInvBPhi3(D, solver, ioms, phi, tol2, max_it2), idx );
  473. } else {
  474. r = vmap(rhs, idx);
  475. }
  476. jsw_timer timer;
  477. errorn = r.norm() / bnrm2;
  478. omega = 1.;
  479. r_tld = r;
  480. Real errornold = 1e14;
  481. // Get down to business
  482. for (iter=0; iter<max_it; ++iter) {
  483. timer.begin();
  484. rho = r_tld.dot(r);
  485. if (abs(rho) < eps) {
  486. ivmap( phi, phi2, idx );
  487. return (0);
  488. }
  489. if (iter > 0) {
  490. beta = (rho/rho_1) * (alpha/omega);
  491. p = r.array() + beta*(p.array()-omega*v.array()).array();
  492. } else {
  493. p = r;
  494. }
  495. tol2 = tol;
  496. max_it2 = 500000;
  497. //v = implicitDCInvBPhi2(D, C, ioms, solver, p, tol2, max_it2);
  498. ivmap(phi, p, idx);
  499. v = vmap(implicitDCInvBPhi3(D, solver, ioms, phi, tol2, max_it2), idx);
  500. alpha = rho / r_tld.dot(v);
  501. s = r.array() - alpha*v.array();
  502. errorn = s.norm()/bnrm2;
  503. if (errorn < tol && iter > 1) {
  504. phi2.array() += alpha*p.array();
  505. ivmap( phi, phi2, idx );
  506. return (0);
  507. }
  508. tol2 = tol;
  509. max_it1 = 500000;
  510. //t = implicitDCInvBPhi2(D, C, ioms, solver, s, tol2, max_it1);
  511. //t = implicitDCInvBPhi3(D, solver, ioms, s, tol2, max_it1);
  512. ivmap(phi, s, idx);
  513. t = vmap(implicitDCInvBPhi3(D, solver, ioms, phi, tol2, max_it1), idx);
  514. omega = t.dot(s) / t.dot(t);
  515. r = s.array() - omega*t.array();
  516. errorn = r.norm() / bnrm2;
  517. iter_done = iter;
  518. if (errorn <= tol) {
  519. ivmap( phi, phi2, idx );
  520. return (0);
  521. }
  522. if (abs(omega) < eps) {
  523. ivmap( phi, phi2, idx );
  524. return (0);
  525. }
  526. rho_1 = rho;
  527. logio << "iteration " << std::setw(3) << iter
  528. << " errorn " << std::setw(6) << std::setprecision(4) << std::scientific << errorn
  529. //<< "\timplicit iterations " << std::setw(5) << max_it1+max_it2
  530. << " time " << std::setw(6) << std::fixed << std::setprecision(2) << timer.end() << std::endl;
  531. // Check to see how progress is going
  532. if (errornold - errorn < 0) {
  533. logio << "Irregular non-monotonic (negative) convergence. Recommend restart. \n";
  534. ivmap( phi, phi2, idx );
  535. return (2);
  536. }
  537. /*
  538. if (errornold - errorn < 1e-14) {
  539. logio << "not making any progress. Giving up\n";
  540. return (1);
  541. }
  542. */
  543. //std::cout << "|| p-s ||" << (alpha*p - omega*s).norm() << std::endl;
  544. // only update phi if good things are happening
  545. phi2.array() += alpha*p.array() + omega*s.array();
  546. errornold = errorn;
  547. }
  548. ivmap( phi, phi2, idx );
  549. return (0);
  550. }
  551. // On Input
  552. // A = Matrix
  553. // B = Right hand side
  554. // X = initial guess, and solution
  555. // maxit = maximum Number of iterations
  556. // tol = error tolerance
  557. // On Output
  558. // X real solution vector
  559. // errorn = Real error norm
  560. template < typename Solver >
  561. int implicitbicgstab_ei(const SparseMat& D,
  562. const Ref< VectorXcr const > ioms,
  563. const Ref< VectorXcr const > rhs,
  564. Ref <VectorXcr> phi,
  565. Solver& solver,
  566. int &max_it, const Real &tol, Real &errorn, int &iter_done, ofstream& logio) {
  567. logio << "using the preconditioned Eigen implicit solver" << std::endl;
  568. Complex omega, rho, rho_1, alpha, beta;
  569. Real tol2;
  570. int iter, max_it2,max_it1;
  571. // Determine size of system and init vectors
  572. int n = phi.size();
  573. VectorXcr r(n);
  574. VectorXcr r_tld(n);
  575. VectorXcr p(n);
  576. VectorXcr v(n);
  577. VectorXcr s(n);
  578. VectorXcr t(n);
  579. // Initialise
  580. iter_done = 0;
  581. Real eps = 1e-100;
  582. Real bnrm2 = rhs.norm();
  583. if (bnrm2 == 0) {
  584. phi.setConstant(0.0);
  585. errorn = 0;
  586. std::cerr << "Trivial case of Ax = b, where b is 0\n";
  587. return (0);
  588. }
  589. // If there is an initial guess
  590. if ( phi.norm() ) {
  591. tol2 = tol;
  592. max_it2 = 50000;
  593. r = rhs - implicitDCInvBPhi3(D, solver, ioms, phi, tol2, max_it2);
  594. } else {
  595. r = rhs;
  596. }
  597. jsw_timer timer;
  598. errorn = r.norm() / bnrm2;
  599. omega = 1.;
  600. r_tld = r;
  601. Real errornold = 1e14;
  602. // Get down to business
  603. for (iter=0; iter<max_it; ++iter) {
  604. timer.begin();
  605. rho = r_tld.dot(r);
  606. if (abs(rho) < eps) return (0);
  607. if (iter > 0) {
  608. beta = (rho/rho_1) * (alpha/omega);
  609. p = r.array() + beta*(p.array()-omega*v.array()).array();
  610. } else {
  611. p = r;
  612. }
  613. tol2 = tol;
  614. max_it2 = 500000;
  615. v = implicitDCInvBPhi3(D, solver, ioms, p, tol2, max_it2);
  616. max_it2 = 0; // solver.iterations();
  617. alpha = rho / r_tld.dot(v);
  618. s = r.array() - alpha*v.array();
  619. errorn = s.norm()/bnrm2;
  620. if (errorn < tol && iter > 1) {
  621. phi.array() += alpha*p.array();
  622. return (0);
  623. }
  624. tol2 = tol;
  625. max_it1 = 500000;
  626. t = implicitDCInvBPhi3(D, solver, ioms, s, tol2, max_it1);
  627. max_it1 = 0; //solver.iterations();
  628. omega = t.dot(s) / t.dot(t);
  629. r = s.array() - omega*t.array();
  630. errorn = r.norm() / bnrm2;
  631. iter_done = iter;
  632. if (errorn <= tol ) return (0);
  633. if (abs(omega) < eps) return (0);
  634. rho_1 = rho;
  635. logio << "iteration " << std::setw(4) << iter
  636. << "\terrorn " << std::setw(6) << std::setprecision(4) << std::scientific << errorn
  637. << "\timplicit iterations " << std::setw(5) << max_it1+max_it2
  638. << "\ttime " << std::setw(10) << std::fixed << std::setprecision(2) << timer.end() << std::endl;
  639. // Check to see how progress is going
  640. if (errornold - errorn < 0) {
  641. logio << "irregular (negative) convergence. Try again? \n";
  642. return (2);
  643. }
  644. // only update phi if good things are happening
  645. phi.array() += alpha*p.array() + omega*s.array();
  646. errornold = errorn;
  647. }
  648. return (0);
  649. }
  650. // On Input
  651. // A = Matrix
  652. // B = Right hand side
  653. // X = initial guess, and solution
  654. // maxit = maximum Number of iterations
  655. // tol = error tolerance
  656. // On Output
  657. // X real solution vector
  658. // errorn = Real error norm
  659. int implicitbicgstab(const SparseMat& D,
  660. const SparseMat& C,
  661. const VectorXcr& ioms,
  662. const SparseMat& MC,
  663. Eigen::Ref< VectorXcr > rhs,
  664. VectorXcr& phi,
  665. int &max_it, Real &tol, Real &errorn, int &iter_done) {
  666. Complex omega, rho, rho_1, alpha, beta;
  667. Real errmin, tol2;
  668. int iter, max_it2;
  669. // // Cholesky decomp
  670. // SparseLLT<SparseMatrix<Complex>, Cholmod>
  671. // CholC(SparseMatrix<Complex> (C.real()) );
  672. // if(!CholC.succeeded()) {
  673. // std::cerr << "decomposiiton failed\n";
  674. // return EXIT_FAILURE;
  675. // }
  676. // Determine size of system and init vectors
  677. int n = phi.size();
  678. VectorXcr r(n);
  679. VectorXcr r_tld(n);
  680. VectorXcr p(n);
  681. VectorXcr v(n);
  682. //VectorXcr p_hat(n);
  683. VectorXcr s(n);
  684. //VectorXcr s_hat(n);
  685. VectorXcr t(n);
  686. VectorXcr xmin(n);
  687. // TODO, refigure for implicit large system
  688. // std::cout << "Start BiCGStab, memory needed: "
  689. // << (sizeof(Complex)*(9+2)*n/(1024.*1024*1024)) << " [Gb]\n";
  690. // Initialise
  691. iter_done = 0;
  692. v.setConstant(0.); // not necessary I don't think
  693. t.setConstant(0.);
  694. Real eps = 1e-100;
  695. Real bnrm2 = rhs.norm();
  696. if (bnrm2 == 0) {
  697. phi.setConstant(0.0);
  698. errorn = 0;
  699. std::cerr << "Trivial case of Ax = b, where b is 0\n";
  700. return (0);
  701. }
  702. // If there is an initial guess
  703. if ( phi.norm() ) {
  704. //r = rhs - A*phi;
  705. tol2 = tol;
  706. max_it2 = 50000;
  707. std::cout << "Initial guess " << std::endl;
  708. r = rhs - implicitDCInvBPhi(D, C, ioms, MC, phi, tol2, max_it2);
  709. //r = rhs - implicitDCInvBPhi (D, C, B, CholC, phi, tol2, max_it2);
  710. } else {
  711. r = rhs;
  712. }
  713. errorn = r.norm() / bnrm2;
  714. //std::cout << "Initial |r| " << r.norm() << "\t" << errorn<< std::endl;
  715. omega = 1.;
  716. r_tld = r;
  717. errmin = 1e30;
  718. Real errornold = 1e6;
  719. // Get down to business
  720. for (iter=0; iter<max_it; ++iter) {
  721. rho = r_tld.dot(r);
  722. if (abs(rho) < eps) return (0);
  723. if (iter > 0) {
  724. beta = (rho/rho_1) * (alpha/omega);
  725. p = r.array() + beta*(p.array()-omega*v.array()).array();
  726. } else {
  727. p = r;
  728. }
  729. // Use pseudo inverse to get approximate answer
  730. //p_hat = p;
  731. tol2 = std::max(1e-4*errorn, tol);
  732. tol2 = tol;
  733. max_it2 = 500000;
  734. //v = A*p_hat;
  735. v = implicitDCInvBPhi(D, C, ioms, MC, p, tol2, max_it2);
  736. //v = implicitDCInvBPhi(D, C, B, CholC, p, tol2, max_it2);
  737. alpha = rho / r_tld.dot(v);
  738. s = r.array() - alpha*v.array();
  739. errorn = s.norm()/bnrm2;
  740. if (errorn < tol && iter > 1) {
  741. phi.array() += alpha*p.array();
  742. return (0);
  743. }
  744. // s_hat = M*s;
  745. //tol2 = tol;
  746. tol2 = std::max(1e-4*errorn, tol);
  747. tol2 = tol;
  748. max_it2 = 50000;
  749. // t = A*s_hat;
  750. t = implicitDCInvBPhi(D, C, ioms, MC, s, tol2, max_it2);
  751. //t = implicitDCInvBPhi(D, C, B, CholC, s, tol2, max_it2);
  752. omega = t.dot(s) / t.dot(t);
  753. phi.array() += alpha*p.array() + omega*s.array();
  754. r = s.array() - omega*t.array();
  755. errorn = r.norm() / bnrm2;
  756. iter_done = iter;
  757. if (errorn < errmin) {
  758. // remember the model with the smallest norm
  759. errmin = errorn;
  760. xmin = phi;
  761. }
  762. if (errorn <= tol ) return (0);
  763. if (abs(omega) < eps) return (0);
  764. rho_1 = rho;
  765. std::cout << "iteration " << std::setw(4) << iter << "\terrorn " << std::setw(6) << std::scientific << errorn
  766. << "\timplicit iterations " << std::setw(5) << max_it2 << std::endl;
  767. if (errornold - errorn < 1e-14) {
  768. std::cout << "not making any progress. Giving up\n";
  769. return (2);
  770. }
  771. errornold = errorn;
  772. }
  773. return (0);
  774. }
  775. //int bicgstab(const SparseMat &A, Eigen::SparseLU< Eigen::SparseMatrix<Complex, RowMajor> ,