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.

EMSchur3D.h 25KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616
  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 02/19/2015 01:10:39 PM
  11. * @version $Id$
  12. * @author Trevor Irons (ti)
  13. * @email Trevor.Irons@xri-geo.com
  14. * @copyright Copyright (c) 2015, XRI Geophysics, LLC
  15. * @copyright Copyright (c) 2015, Trevor Irons
  16. * @copyright Copyright (c) 2011, Trevor Irons
  17. * @copyright Copyright (c) 2011, Colorado School of Mines
  18. */
  19. #ifndef EMSCHUR3D_INC
  20. #define EMSCHUR3D_INC
  21. #include "EMSchur3DBase.h"
  22. #include "CSymSimplicialCholesky.h"
  23. namespace Lemma {
  24. /**
  25. \brief Templated concrete classes of EMSChur3DBase.
  26. \details
  27. */
  28. template < class Solver >
  29. class EMSchur3D : public EMSchur3DBase {
  30. //friend std::ostream &operator<<(std::ostream &stream,
  31. // const EMSchur3D &ob);
  32. public:
  33. // ==================== LIFECYCLE =======================
  34. /**
  35. * @copybrief LemmaObject::New()
  36. * @copydetails LemmaObject::New()
  37. */
  38. static EMSchur3D* New() {
  39. EMSchur3D<Solver>* Obj = new EMSchur3D<Solver>("EMSchur3D");
  40. Obj->AttachTo(Obj);
  41. return Obj;
  42. }
  43. /**
  44. * @copybrief LemmaObject::Delete()
  45. * @copydetails LemmaObject::Delete()
  46. */
  47. void Delete() {
  48. this->DetachFrom(this);
  49. }
  50. // ==================== OPERATORS =======================
  51. // ==================== OPERATIONS =======================
  52. /** Solves a single source problem. This method is thread safe.
  53. * @param[in] Source is the source term for generating primary fields
  54. * @param[in] isource is the source index
  55. */
  56. void SolveSource( DipoleSource* Source , const int& isource);
  57. /** Builds the solver for the C matrix */
  58. void BuildCDirectSolver( );
  59. // ==================== ACCESS =======================
  60. // ==================== INQUIRY =======================
  61. #ifdef HAVE_YAMLCPP
  62. // /**
  63. // * Uses YAML to serialize this object.
  64. // * @return a YAML::Node
  65. // */
  66. // YAML::Node Serialize() const;
  67. //
  68. // /**
  69. // * Constructs an object from a YAML::Node.
  70. // */
  71. // static EMSchur3D* DeSerialize(const YAML::Node& node);
  72. #endif
  73. protected:
  74. // ==================== LIFECYCLE =======================
  75. /** Default protected constructor, use New */
  76. EMSchur3D (const std::string& name) : EMSchur3DBase(name), CSolver(NULL) {
  77. }
  78. // #ifdef HAVE_YAMLCPP
  79. // /** Protected DeDerializing constructor, use factory DeSerialize method*/
  80. // EMSchur3D (const YAML::Node& node): EMSchur3DBase(node), CSolver(NULL) {
  81. // }
  82. // #endif
  83. /** Default protected destructor, use Delete */
  84. ~EMSchur3D () {
  85. // TODO delete arrays
  86. }
  87. /**
  88. * @copybrief LemmaObject::Release()
  89. * @copydetails LemmaObject::Release()
  90. */
  91. void Release() {
  92. delete this;
  93. }
  94. private:
  95. // ==================== DATA MEMBERS =========================
  96. /** The templated solver for C */
  97. Solver* CSolver;
  98. Eigen::SparseMatrix<Complex> Csym;
  99. }; // ----- end of class EMSchur3D -----
  100. ////////////////////////////////////////////////////////////////////////////////////////
  101. // Implimentation and Specialisations //
  102. ////////////////////////////////////////////////////////////////////////////////////////
  103. //--------------------------------------------------------------------------------------
  104. // Class: EMSchur3D
  105. // Method: SolveSource
  106. //--------------------------------------------------------------------------------------
  107. template < class Solver >
  108. void EMSchur3D<Solver>::SolveSource ( DipoleSource* Source, const int& isource ) {
  109. // figure out which omega we are working with
  110. int iw = -1;
  111. for (int iiw=0; iiw<Omegas.size(); ++iiw) {
  112. if (Omegas[iiw] - Source->GetAngularFrequency(0) < 1e-3 ) {
  113. iw = iiw;
  114. }
  115. }
  116. if (iw == -1) {
  117. std::cerr << "FREQUENCY DOOM IN EMSchur3D::SolveSource \n";
  118. exit(EXIT_FAILURE);
  119. }
  120. ///////////////////////////////////
  121. // Set up primary fields
  122. // TODO, this is a little stupid as they all share the same points. We need to extend
  123. // EmEARTH to be able to input a grid so that points are not explicitly needed like
  124. // this. This requires some care as calcs are made on faces.
  125. // Alternatively, the bins function of ReceiverPoints could be extended quite easily.
  126. // This may be the way to do this.
  127. Lemma::ReceiverPoints* dpoint = Lemma::ReceiverPoints::New();
  128. FillPoints(dpoint);
  129. PrimaryField(Source, dpoint);
  130. // Allocate a ton of memory
  131. VectorXcr Phi = VectorXcr::Zero(uns);
  132. VectorXcr ms(unx+uny+unz); // mu sigma
  133. // Vector potential (A) Vector and phi
  134. VectorXcr Se = VectorXcr::Zero(unx+uny+unz);
  135. //VectorXcr A = VectorXcr::Zero(unx+uny+unz);
  136. VectorXcr E = VectorXcr::Zero(unx+uny+unz);
  137. VectorXcr E0 = VectorXcr::Zero(unx+uny+unz);
  138. // Lets get cracking
  139. FillSourceTerms(ms, Se, E0, dpoint, Omegas[iw]);
  140. /////////////////////////////////////////////////
  141. // LOG File
  142. std::string logfile (ResFile);
  143. logfile += to_string(isource) + std::string(".log");
  144. ofstream logio(logfile.c_str());
  145. logio << *Source << std::endl;
  146. logio << *Grid << std::endl;
  147. logio << *LayModel << std::endl;
  148. // solve for RHS
  149. int max_it(nx*ny*nz), iter_done(0);
  150. Real tol(3e-16), errorn(0);
  151. logio << "solving RHS for source " << isource << std::endl;
  152. // TODO, this is stupid, try and get rid of this copy!
  153. Eigen::SparseMatrix<Complex> Cc = Cvec[iw];
  154. jsw_timer timer;
  155. jsw_timer timer2;
  156. timer.begin();
  157. timer2.begin();
  158. /////////////////////////////////////////
  159. // Solve for RHS
  160. VectorXcr A = CSolver[iw].solve(Se);
  161. // // Solve Real system instead
  162. // The Real system is quasi-definite, though an LDLT decomposition exists, CHOLMOD doesn't find it.
  163. // An LU can be done on this, but compute performance is very similiar to the complex system, and diagonal pivoting
  164. // cannot be assumed to be best, hurting solve time.
  165. // /* EXPERIMENTAL */
  166. // VectorXr b2 = VectorXr::Zero(2*(unx+uny+unz));
  167. // b2.head(unx+uny+unz) = Se.real();
  168. // b2.tail(unx+uny+unz) = Se.imag();
  169. // VectorXr A2 = CReSolver[iw].solve(b2);
  170. // A.real() = A2.head( unx+uny+unz );
  171. // A.imag() = -A2.tail( unx+uny+unz ); // Due to decomp. negative!
  172. // /* END EXPERIMENTAL */
  173. VectorXcr ADiv = D*A; // ADiv == RHS == D C^I Se
  174. VectorXcr Error = ((Cc.selfadjointView<Eigen::Lower>()*A).array() - Se.array());
  175. logio << "|| Div(A) || = " << ADiv.norm()
  176. // << " in " << iter_done << " iterations"
  177. //<< " with error " << errorn << "\t"
  178. << "\tInital solution error "<< Error.norm() // Iteritive info
  179. << "\ttime " << timer.end() << std::endl;
  180. //VectorXcr ADivMAC = ADiv.array() * MAC.array().cast<Complex>();
  181. //logio << "|| Div(A) || on MAC grid " << ADivMAC.norm() << std::endl;
  182. /////////////////////
  183. // Solve for Phi
  184. logio << "Solving for Phi " << std::flush;
  185. timer.begin();
  186. tol = 1e-18;
  187. int success(2);
  188. success = implicitbicgstab(D, idx, ms, ADiv, Phi, CSolver[iw], max_it, tol, errorn, iter_done, logio);
  189. //Phi.array() *= MAC.array().cast<Complex>(); // remove phi from air regions
  190. /* Restart if necessary */
  191. int nrestart(1);
  192. // TODO send MAC to implicitbicgstab?
  193. while (success == 2 && nrestart < 18 && iter_done > 1) {
  194. success = implicitbicgstab(D, idx, ms, ADiv, Phi, CSolver[iw], max_it, tol, errorn, iter_done, logio);
  195. //Phi.array() *= MAC.array().cast<Complex>(); // remove phi from air regions
  196. nrestart += 1;
  197. }
  198. logio << "Implicit BiCGStab solution in " << iter_done << " iterations."
  199. << " with error " << std::setprecision(8) << std::scientific << errorn << std::endl;
  200. logio << "time "<< timer.end() << " [s]" << std::endl;
  201. E = ms.array()*(D.transpose()*Phi).array(); // Temp, field due to charge
  202. /////////////////////////////////////
  203. // Compute A
  204. /////////////////////////////////////
  205. logio << "Solving for A using phi" << std::endl;
  206. std::cout << "Solving for A" << std::endl;
  207. max_it = nx*ny*nz;
  208. tol = 5e-16;
  209. errorn = 0;
  210. iter_done = 0;
  211. timer.begin();
  212. A = CSolver[iw].solve( (Se-E).eval() ); // UmfPack requires eval?
  213. VectorXcr ADiv2 = D*A;
  214. logio << "|| Div(A) || = " << ADiv2.norm() ;
  215. //" in " << iter_done << " iterations"
  216. //<< " with error " << errorn << "\t";
  217. // Report error of solutions
  218. Error = ((Cc.selfadjointView<Eigen::Lower>()*A).array() + E.array() - Se.array());
  219. logio << "\tsolution error " << Error.norm()
  220. << std::fixed << std::setprecision(2) << "\ttime " << timer.end() << "\ttotal time " << timer2.end() << std::endl;
  221. logio.close();
  222. //////////////////////////////////////
  223. // Update Fields and report
  224. E.array() = Complex(0,-Omegas[iw])*A.array() - (D.transpose()*Phi).array(); // Secondary Field Only
  225. VectorXcr B = StaggeredGridCurl(A);
  226. WriteVTKResults( ResFile+ to_string(isource), A, Se, E0, E , Phi, ADiv, ADiv2, B);
  227. dpoint->Delete();
  228. return ;
  229. } // ----- end of method EMSchur3D::SolveSource -----
  230. //--------------------------------------------------------------------------------------
  231. // Class: EMSchur3DBase
  232. // Method: BuildCDirectSolver
  233. //--------------------------------------------------------------------------------------
  234. template < class Solver >
  235. void EMSchur3D<Solver>::BuildCDirectSolver ( ) {
  236. CSolver = new Solver[Omegas.size()];
  237. for (int iw=0; iw<Omegas.size(); ++iw) {
  238. jsw_timer timer;
  239. timer.begin();
  240. /* Complex system */
  241. /*
  242. std::cout << "Generic solver pattern analyzing C_" << iw << ",";
  243. std::cout.flush();
  244. CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  245. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  246. // factorize
  247. timer.begin();
  248. std::cout << "Generic solver factorising C_" << iw << ", ";
  249. std::cout.flush();
  250. CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  251. */
  252. std::cerr << "No solver Specified!" << iw << ",";
  253. exit(EXIT_FAILURE);
  254. //CSolver[iw].compute( Cvec[iw].selfadjointView< Eigen::Lower>() );
  255. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  256. }
  257. }
  258. #ifdef HAVE_SUPERLUMT
  259. template<>
  260. void EMSchur3D< Eigen::SuperLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::BuildCDirectSolver() {
  261. CSolver = new Eigen::SuperLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > [Omegas.size()];
  262. for (int iw=0; iw<Omegas.size(); ++iw) {
  263. jsw_timer timer;
  264. timer.begin();
  265. /* SuperLU */
  266. //CSolver[iw].options().DiagPivotThresh = 0.01;
  267. //CSolver[iw].options().SymmetricMode = YES;
  268. //CSolver[iw].options().ColPerm = MMD_AT_PLUS_A;
  269. //CSolver[iw].options().Trans = NOTRANS;
  270. //CSolver[iw].options().ConditionNumber = NO;
  271. //std::cout << "SuperLU options:\n";
  272. //std::cout << "\tPivot Threshold: " << CSolver[iw].options().DiagPivotThresh << std::endl;
  273. //std::cout << "\tSymmetric mode: " << CSolver[iw].options().SymmetricMode << std::endl;
  274. //std::cout << "\tEquilibrate: " << CSolver[iw].options().Equil << std::endl;
  275. //std::cout << "\tCol Permutation: " << CSolver[iw].options().ColPerm << std::endl;
  276. //std::cout << "\tTrans: " << CSolver[iw].options().Trans << std::endl;
  277. //std::cout << "\tCondition Number: " << CSolver[iw].options().ConditionNumber << std::endl;
  278. /* Complex system */
  279. std::cout << "SuperLU_MT pattern analyzing C_" << iw << ",";
  280. std::cout.flush();
  281. CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  282. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  283. // factorize
  284. timer.begin();
  285. std::cout << "SuperLU_MT factorising C_" << iw << ", ";
  286. std::cout.flush();
  287. CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  288. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  289. }
  290. }
  291. #endif
  292. template<>
  293. void EMSchur3D< Eigen::SparseLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::COLAMDOrdering<int> > >::BuildCDirectSolver() {
  294. CSolver = new Eigen::SparseLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::COLAMDOrdering<int> > [Omegas.size()];
  295. for (int iw=0; iw<Omegas.size(); ++iw) {
  296. jsw_timer timer;
  297. timer.begin();
  298. CSolver[iw].isSymmetric(true);
  299. CSolver[iw].setPivotThreshold(0.0);
  300. /* Complex system */
  301. std::cout << "SparseLU pattern analyzing C_" << iw << ",";
  302. std::cout.flush();
  303. CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  304. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  305. // factorize
  306. timer.begin();
  307. std::cout << "SparseLU factorising C_" << iw << ", ";
  308. std::cout.flush();
  309. CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  310. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  311. }
  312. }
  313. // template<>
  314. // void EMSchur3D< Eigen::CholmodSupernodalLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
  315. // CSolver = new Eigen::CholmodSupernodalLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > [Omegas.size()];
  316. // for (int iw=0; iw<Omegas.size(); ++iw) {
  317. // Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  318. // jsw_timer timer;
  319. // timer.begin();
  320. // /* Complex system */
  321. // std::cout << "CholmodSupernodalLLT pattern analyzing C_" << iw << ",";
  322. // std::cout.flush();
  323. // CSolver[iw].analyzePattern( Csym );
  324. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  325. // /* factorize */
  326. // timer.begin();
  327. // std::cout << "CholmodSupernodalLLT factorising C_" << iw << ", ";
  328. // std::cout.flush();
  329. // CSolver[iw].factorize( Csym );
  330. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  331. // }
  332. // }
  333. template<>
  334. void EMSchur3D< Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::NaturalOrdering<int> > > ::BuildCDirectSolver() {
  335. CSolver = new Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::NaturalOrdering<int> > [Omegas.size()];
  336. for (int iw=0; iw<Omegas.size(); ++iw) {
  337. Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  338. jsw_timer timer;
  339. timer.begin();
  340. /* Complex system */
  341. std::cout << "CSymSimplicialLLT<NaturalOrdering> pattern analyzing C_" << iw << ",";
  342. std::cout.flush();
  343. CSolver[iw].analyzePattern( Csym );
  344. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  345. /* factorize */
  346. timer.begin();
  347. std::cout << "CSymSimplicialLLT<NaturalOrdering> factorising C_" << iw << ", ";
  348. std::cout.flush();
  349. CSolver[iw].factorize( Csym );
  350. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  351. }
  352. }
  353. template<>
  354. void EMSchur3D< Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > > ::BuildCDirectSolver() {
  355. CSolver = new Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > [Omegas.size()];
  356. for (int iw=0; iw<Omegas.size(); ++iw) {
  357. //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  358. jsw_timer timer;
  359. timer.begin();
  360. /* Complex system */
  361. std::cout << "CSymSimplicialLLT<AMDOrdering> pattern analyzing C_" << iw << ",";
  362. std::cout.flush();
  363. CSolver[iw].analyzePattern( Cvec[iw] );
  364. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  365. /* factorize */
  366. timer.begin();
  367. std::cout << "CSymSimplicialLLT<AMDOrdering> factorising C_" << iw << ", ";
  368. std::cout.flush();
  369. CSolver[iw].factorize( Cvec[iw] );
  370. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  371. }
  372. }
  373. template<>
  374. void EMSchur3D< Eigen::CSymSimplicialLDLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > > ::BuildCDirectSolver() {
  375. CSolver = new Eigen::CSymSimplicialLDLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > [Omegas.size()];
  376. for (int iw=0; iw<Omegas.size(); ++iw) {
  377. Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  378. jsw_timer timer;
  379. timer.begin();
  380. /* Complex system */
  381. std::cout << "CSymSimplicialLDLT<AMDOrdering> pattern analyzing C_" << iw << ",";
  382. std::cout.flush();
  383. CSolver[iw].analyzePattern( Csym );
  384. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  385. /* factorize */
  386. timer.begin();
  387. std::cout << "CSymSimplicialLDLT<AMDOrdering> factorising C_" << iw << ", ";
  388. std::cout.flush();
  389. CSolver[iw].factorize( Csym );
  390. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  391. }
  392. }
  393. template<>
  394. void EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::IncompleteLUT<Complex> > > ::BuildCDirectSolver() {
  395. CSolver = new Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::IncompleteLUT<Complex> > [Omegas.size()];
  396. for (int iw=0; iw<Omegas.size(); ++iw) {
  397. Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  398. jsw_timer timer;
  399. timer.begin();
  400. /* Complex system */
  401. std::cout << "BiCGSTAB(ILU) pattern analyzing C_" << iw << ",";
  402. std::cout.flush();
  403. CSolver[iw].analyzePattern( Csym );
  404. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  405. /* factorize */
  406. timer.begin();
  407. std::cout << "BiCGSTAB(ILU) factorising C_" << iw << ", ";
  408. std::cout.flush();
  409. CSolver[iw].factorize( Csym );
  410. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  411. }
  412. }
  413. template<>
  414. void EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > > ::BuildCDirectSolver() {
  415. CSolver = new Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > [Omegas.size()];
  416. for (int iw=0; iw<Omegas.size(); ++iw) {
  417. Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  418. jsw_timer timer;
  419. timer.begin();
  420. /* Complex system */
  421. std::cout << "BiCGSTAB pattern analyzing C_" << iw << ",";
  422. std::cout.flush();
  423. CSolver[iw].analyzePattern( Csym );
  424. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  425. // factorize
  426. timer.begin();
  427. std::cout << "BiCGSTAB factorising C_" << iw << ", ";
  428. std::cout.flush();
  429. CSolver[iw].factorize( Csym );
  430. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  431. }
  432. }
  433. template<>
  434. void EMSchur3D< Eigen::ConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
  435. CSolver = new Eigen::ConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > [Omegas.size()];
  436. for (int iw=0; iw<Omegas.size(); ++iw) {
  437. //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  438. jsw_timer timer;
  439. timer.begin();
  440. /* Complex system */
  441. std::cout << "ConjugateGradient pattern analyzing C_" << iw << ",";
  442. std::cout.flush();
  443. CSolver[iw].analyzePattern( Cvec[iw] );
  444. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  445. // factorize
  446. timer.begin();
  447. std::cout << "ConjugateGradient factorising C_" << iw << ", ";
  448. std::cout.flush();
  449. CSolver[iw].factorize( Cvec[iw] );
  450. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  451. }
  452. }
  453. // template<>
  454. // void EMSchur3D< Eigen::PastixLLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
  455. // CSolver = new Eigen::PastixLLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > [Omegas.size()];
  456. // //MPI_Init(NULL, NULL);
  457. // for (int iw=0; iw<Omegas.size(); ++iw) {
  458. // //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  459. // jsw_timer timer;
  460. // timer.begin();
  461. // /* Complex system */
  462. // std::cout << "PaStiX LLT pattern analyzing C_" << iw << ",";
  463. // std::cout.flush();
  464. // CSolver[iw].analyzePattern( Cvec[iw] );
  465. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  466. // // factorize
  467. // timer.begin();
  468. // std::cout << "PaStiX LLT factorising C_" << iw << ", ";
  469. // std::cout.flush();
  470. // CSolver[iw].factorize( Cvec[iw] );
  471. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  472. // }
  473. // }
  474. //
  475. // template<>
  476. // void EMSchur3D< Eigen::PastixLDLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
  477. // CSolver = new Eigen::PastixLDLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > [Omegas.size()];
  478. // //MPI_Init(NULL, NULL);
  479. // for (int iw=0; iw<Omegas.size(); ++iw) {
  480. // //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  481. // jsw_timer timer;
  482. // timer.begin();
  483. // /* Complex system */
  484. // std::cout << "PaStiX LDLT pattern analyzing C_" << iw << ",";
  485. // std::cout.flush();
  486. // CSolver[iw].analyzePattern( Cvec[iw] );
  487. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  488. // // factorize
  489. // timer.begin();
  490. // std::cout << "PaStiX LDLT factorising C_" << iw << ", ";
  491. // std::cout.flush();
  492. // CSolver[iw].factorize( Cvec[iw] );
  493. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  494. // std::cout << "INFO " << CSolver[iw].info( ) << std::endl;
  495. // }
  496. // }
  497. //
  498. // template<>
  499. // void EMSchur3D< Eigen::PastixLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, true > > ::BuildCDirectSolver() {
  500. // CSolver = new Eigen::PastixLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, true > [Omegas.size()];
  501. // //MPI_Init(NULL, NULL);
  502. // for (int iw=0; iw<Omegas.size(); ++iw) {
  503. // Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  504. // jsw_timer timer;
  505. // timer.begin();
  506. // /* Complex system */
  507. // std::cout << "PaStiX LU pattern analyzing C_" << iw << ",";
  508. // std::cout.flush();
  509. // CSolver[iw].compute( Csym );
  510. // std::cout << "PaStiX LU Done C_" << iw << std::endl;;
  511. // // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  512. // // // factorize
  513. // // timer.begin();
  514. // // std::cout << "PaStiX LU factorising C_" << iw << ", ";
  515. // // std::cout.flush();
  516. // // CSolver[iw].factorize( Csym );
  517. // // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  518. // }
  519. // }
  520. } // ----- end of Lemma name -----
  521. #endif // ----- #ifndef EMSCHUR3D_INC -----