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 26KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622
  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 "bicgstab.h"
  23. //#include "CSymSimplicialCholesky.h"
  24. namespace Lemma {
  25. /**
  26. \brief Templated concrete classes of EMSChur3DBase.
  27. \details
  28. */
  29. template < class Solver >
  30. class EMSchur3D : public EMSchur3DBase {
  31. friend std::ostream &operator << (std::ostream &stream, const EMSchur3D &ob) {
  32. stream << ob.Serialize() << "\n---\n"; // End of doc
  33. return stream;
  34. }
  35. struct ctor_key {};
  36. //friend std::ostream &operator<<(std::ostream &stream,
  37. // const EMSchur3D &ob);
  38. public:
  39. // ==================== LIFECYCLE =======================
  40. /**
  41. * @copybrief LemmaObject::New()
  42. * @copydetails LemmaObject::New()
  43. */
  44. static std::shared_ptr< EMSchur3D > NewSP() {
  45. return std::make_shared< EMSchur3D<Solver> >( ctor_key() );
  46. //return std::make_shared< EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > > >( ctor_key() ) ;
  47. }
  48. /** Default protected constructor, use New */
  49. explicit EMSchur3D ( const ctor_key& key ) : EMSchur3DBase( ), CSolver( nullptr ) {
  50. }
  51. /** Locked DeDerializing constructor, use factory DeSerialize method*/
  52. EMSchur3D (const YAML::Node& node): EMSchur3DBase(node), CSolver( nullptr ) {
  53. }
  54. /** Default protected destructor, use Delete */
  55. virtual ~EMSchur3D () {
  56. // TODO delete arrays
  57. }
  58. /**
  59. * Uses YAML to serialize this object.
  60. * @return a YAML::Node
  61. */
  62. YAML::Node Serialize() const {
  63. YAML::Node node = EMSchur3DBase::Serialize();
  64. //node["NumberOfLayers"] = NumberOfLayers;
  65. node.SetTag( this->GetName() );
  66. return node;
  67. }
  68. /**
  69. * Constructs an object from a YAML::Node.
  70. */
  71. static EMSchur3D* DeSerialize(const YAML::Node& node);
  72. // ==================== OPERATORS =======================
  73. // ==================== OPERATIONS =======================
  74. /** Solves a single source problem. This method is thread safe.
  75. * @param[in] Source is the source term for generating primary fields
  76. * @param[in] isource is the source index
  77. */
  78. void SolveSource( std::shared_ptr<DipoleSource> Source , const int& isource);
  79. /** Builds the solver for the C matrix */
  80. void BuildCDirectSolver( );
  81. // ==================== ACCESS =======================
  82. virtual std::string GetName() const {
  83. return this->CName;
  84. }
  85. // ==================== INQUIRY =======================
  86. protected:
  87. // ==================== LIFECYCLE =======================
  88. private:
  89. /** Copy constructor */
  90. EMSchur3D( const EMSchur3D& ) = delete;
  91. // ==================== DATA MEMBERS =========================
  92. /** The templated solver for C */
  93. Solver* CSolver;
  94. Eigen::SparseMatrix<Complex> Csym;
  95. static constexpr auto CName = "EMSchur3D";
  96. }; // ----- end of class EMSchur3D -----
  97. ////////////////////////////////////////////////////////////////////////////////////////
  98. // Implimentation and Specialisations //
  99. ////////////////////////////////////////////////////////////////////////////////////////
  100. //--------------------------------------------------------------------------------------
  101. // Class: EMSchur3D
  102. // Method: SolveSource
  103. //--------------------------------------------------------------------------------------
  104. template < class Solver >
  105. void EMSchur3D<Solver>::SolveSource ( std::shared_ptr<DipoleSource> Source, const int& isource ) {
  106. // figure out which omega we are working with
  107. int iw = -1;
  108. for (int iiw=0; iiw<Omegas.size(); ++iiw) {
  109. if (Omegas[iiw] - Source->GetAngularFrequency(0) < 1e-3 ) {
  110. iw = iiw;
  111. }
  112. }
  113. if (iw == -1) {
  114. std::cerr << "FREQUENCY DOOM IN EMSchur3D::SolveSource \n";
  115. exit(EXIT_FAILURE);
  116. }
  117. ///////////////////////////////////
  118. // Set up primary fields
  119. // TODO, this is a little stupid as they all share the same points. We need to extend
  120. // EmEARTH to be able to input a grid so that points are not explicitly needed like
  121. // this. This requires some care as calcs are made on faces.
  122. // Alternatively, the bins function of ReceiverPoints could be extended quite easily.
  123. // This may be the way to do this.
  124. //Lemma::ReceiverPoints* dpoint = Lemma::ReceiverPoints::New();
  125. std::shared_ptr< FieldPoints > dpoint = FieldPoints::NewSP();
  126. FillPoints(dpoint);
  127. PrimaryField(Source, dpoint);
  128. // Allocate a ton of memory
  129. VectorXcr Phi = VectorXcr::Zero(uns);
  130. VectorXcr ms(unx+uny+unz); // mu sigma
  131. // Vector potential (A) Vector and phi
  132. VectorXcr Se = VectorXcr::Zero(unx+uny+unz);
  133. //VectorXcr A = VectorXcr::Zero(unx+uny+unz);
  134. VectorXcr E = VectorXcr::Zero(unx+uny+unz);
  135. VectorXcr E0 = VectorXcr::Zero(unx+uny+unz);
  136. // Lets get cracking
  137. FillSourceTerms(ms, Se, E0, dpoint, Omegas[iw]);
  138. /////////////////////////////////////////////////
  139. // LOG File
  140. std::string logfile (ResFile);
  141. logfile += to_string(isource) + std::string(".log");
  142. ofstream logio(logfile.c_str());
  143. logio << *Source << std::endl;
  144. logio << *Grid << std::endl;
  145. logio << *LayModel << std::endl;
  146. // solve for RHS
  147. int max_it(nx*ny*nz), iter_done(0);
  148. Real tol(3e-16), errorn(0);
  149. logio << "solving RHS for source " << isource << std::endl;
  150. // TODO, this is stupid, try and get rid of this copy!
  151. Eigen::SparseMatrix<Complex> Cc = Cvec[iw];
  152. jsw_timer timer;
  153. jsw_timer timer2;
  154. timer.begin();
  155. timer2.begin();
  156. /////////////////////////////////////////
  157. // Solve for RHS
  158. VectorXcr A = CSolver[iw].solve(Se);
  159. // // Solve Real system instead
  160. // The Real system is quasi-definite, though an LDLT decomposition exists, CHOLMOD doesn't find it.
  161. // An LU can be done on this, but compute performance is very similiar to the complex system, and diagonal pivoting
  162. // cannot be assumed to be best, hurting solve time.
  163. // /* EXPERIMENTAL */
  164. // VectorXr b2 = VectorXr::Zero(2*(unx+uny+unz));
  165. // b2.head(unx+uny+unz) = Se.real();
  166. // b2.tail(unx+uny+unz) = Se.imag();
  167. // VectorXr A2 = CReSolver[iw].solve(b2);
  168. // A.real() = A2.head( unx+uny+unz );
  169. // A.imag() = -A2.tail( unx+uny+unz ); // Due to decomp. negative!
  170. // /* END EXPERIMENTAL */
  171. VectorXcr ADiv = D*A; // ADiv == RHS == D C^I Se
  172. VectorXcr Error = ((Cc.selfadjointView<Eigen::Lower>()*A).array() - Se.array());
  173. logio << "|| Div(A) || = " << ADiv.norm()
  174. // << " in " << iter_done << " iterations"
  175. //<< " with error " << errorn << "\t"
  176. << "\tInital solution error "<< Error.norm() // Iteritive info
  177. << "\ttime " << timer.end() << std::endl;
  178. //VectorXcr ADivMAC = ADiv.array() * MAC.array().cast<Complex>();
  179. //logio << "|| Div(A) || on MAC grid " << ADivMAC.norm() << std::endl;
  180. /////////////////////
  181. // Solve for Phi
  182. logio << "Solving for Phi " << std::flush;
  183. timer.begin();
  184. tol = 1e-18;
  185. int success(2);
  186. success = implicitbicgstab(D, idx, ms, ADiv, Phi, CSolver[iw], max_it, tol, errorn, iter_done, logio);
  187. //Phi.array() *= MAC.array().cast<Complex>(); // remove phi from air regions
  188. /* Restart if necessary */
  189. int nrestart(1);
  190. // TODO send MAC to implicitbicgstab?
  191. while (success == 2 && nrestart < 18 && iter_done > 1) {
  192. success = implicitbicgstab(D, idx, ms, ADiv, Phi, CSolver[iw], max_it, tol, errorn, iter_done, logio);
  193. //Phi.array() *= MAC.array().cast<Complex>(); // remove phi from air regions
  194. nrestart += 1;
  195. }
  196. logio << "Implicit BiCGStab solution in " << iter_done << " iterations."
  197. << " with error " << std::setprecision(8) << std::scientific << errorn << std::endl;
  198. logio << "time "<< timer.end() << " [s]" << std::endl;
  199. E = ms.array()*(D.transpose()*Phi).array(); // Temp, field due to charge
  200. /////////////////////////////////////
  201. // Compute A
  202. /////////////////////////////////////
  203. logio << "Solving for A using phi" << std::endl;
  204. std::cout << "Solving for A" << std::endl;
  205. max_it = nx*ny*nz;
  206. tol = 5e-16;
  207. errorn = 0;
  208. iter_done = 0;
  209. timer.begin();
  210. A = CSolver[iw].solve( (Se-E).eval() ); // UmfPack requires eval?
  211. VectorXcr ADiv2 = D*A;
  212. logio << "|| Div(A) || = " << ADiv2.norm() ;
  213. //" in " << iter_done << " iterations"
  214. //<< " with error " << errorn << "\t";
  215. // Report error of solutions
  216. Error = ((Cc.selfadjointView<Eigen::Lower>()*A).array() + E.array() - Se.array());
  217. logio << "\tsolution error " << Error.norm()
  218. << std::fixed << std::setprecision(2) << "\ttime " << timer.end() << "\ttotal time " << timer2.end() << std::endl;
  219. logio.close();
  220. //////////////////////////////////////
  221. // Update Fields and report
  222. E.array() = Complex(0,-Omegas[iw])*A.array() - (D.transpose()*Phi).array(); // Secondary Field Only
  223. VectorXcr B = StaggeredGridCurl(A);
  224. WriteVTKResults( ResFile+ to_string(isource), A, Se, E0, E , Phi, ADiv, ADiv2, B);
  225. //dpoint->Delete();
  226. return ;
  227. } // ----- end of method EMSchur3D::SolveSource -----
  228. //--------------------------------------------------------------------------------------
  229. // Class: EMSchur3DBase
  230. // Method: BuildCDirectSolver
  231. //--------------------------------------------------------------------------------------
  232. template < class Solver >
  233. void EMSchur3D<Solver>::BuildCDirectSolver ( ) {
  234. CSolver = new Solver[Omegas.size()];
  235. for (int iw=0; iw<Omegas.size(); ++iw) {
  236. jsw_timer timer;
  237. timer.begin();
  238. /* Complex system */
  239. /*
  240. std::cout << "Generic solver pattern analyzing C_" << iw << ",";
  241. std::cout.flush();
  242. CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  243. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  244. // factorize
  245. timer.begin();
  246. std::cout << "Generic solver factorising C_" << iw << ", ";
  247. std::cout.flush();
  248. CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  249. */
  250. std::cerr << "No solver Specified!" << iw << ",";
  251. exit(EXIT_FAILURE);
  252. //CSolver[iw].compute( Cvec[iw].selfadjointView< Eigen::Lower>() );
  253. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  254. }
  255. }
  256. #ifdef HAVE_SUPERLUMT
  257. template<>
  258. void EMSchur3D< Eigen::SuperLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::BuildCDirectSolver() {
  259. CSolver = new Eigen::SuperLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > [Omegas.size()];
  260. for (int iw=0; iw<Omegas.size(); ++iw) {
  261. jsw_timer timer;
  262. timer.begin();
  263. /* SuperLU */
  264. //CSolver[iw].options().DiagPivotThresh = 0.01;
  265. //CSolver[iw].options().SymmetricMode = YES;
  266. //CSolver[iw].options().ColPerm = MMD_AT_PLUS_A;
  267. //CSolver[iw].options().Trans = NOTRANS;
  268. //CSolver[iw].options().ConditionNumber = NO;
  269. //std::cout << "SuperLU options:\n";
  270. //std::cout << "\tPivot Threshold: " << CSolver[iw].options().DiagPivotThresh << std::endl;
  271. //std::cout << "\tSymmetric mode: " << CSolver[iw].options().SymmetricMode << std::endl;
  272. //std::cout << "\tEquilibrate: " << CSolver[iw].options().Equil << std::endl;
  273. //std::cout << "\tCol Permutation: " << CSolver[iw].options().ColPerm << std::endl;
  274. //std::cout << "\tTrans: " << CSolver[iw].options().Trans << std::endl;
  275. //std::cout << "\tCondition Number: " << CSolver[iw].options().ConditionNumber << std::endl;
  276. /* Complex system */
  277. std::cout << "SuperLU_MT pattern analyzing C_" << iw << ",";
  278. std::cout.flush();
  279. CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  280. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  281. // factorize
  282. timer.begin();
  283. std::cout << "SuperLU_MT factorising C_" << iw << ", ";
  284. std::cout.flush();
  285. CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  286. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  287. }
  288. }
  289. #endif
  290. template<>
  291. void EMSchur3D< Eigen::SparseLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::COLAMDOrdering<int> > >::BuildCDirectSolver() {
  292. CSolver = new Eigen::SparseLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::COLAMDOrdering<int> > [Omegas.size()];
  293. for (int iw=0; iw<Omegas.size(); ++iw) {
  294. jsw_timer timer;
  295. timer.begin();
  296. CSolver[iw].isSymmetric(true);
  297. CSolver[iw].setPivotThreshold(0.0);
  298. /* Complex system */
  299. std::cout << "SparseLU pattern analyzing C_" << iw << ",";
  300. std::cout.flush();
  301. CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  302. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  303. // factorize
  304. timer.begin();
  305. std::cout << "SparseLU factorising C_" << iw << ", ";
  306. std::cout.flush();
  307. CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  308. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  309. }
  310. }
  311. // template<>
  312. // void EMSchur3D< Eigen::CholmodSupernodalLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
  313. // CSolver = new Eigen::CholmodSupernodalLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > [Omegas.size()];
  314. // for (int iw=0; iw<Omegas.size(); ++iw) {
  315. // Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  316. // jsw_timer timer;
  317. // timer.begin();
  318. // /* Complex system */
  319. // std::cout << "CholmodSupernodalLLT pattern analyzing C_" << iw << ",";
  320. // std::cout.flush();
  321. // CSolver[iw].analyzePattern( Csym );
  322. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  323. // /* factorize */
  324. // timer.begin();
  325. // std::cout << "CholmodSupernodalLLT factorising C_" << iw << ", ";
  326. // std::cout.flush();
  327. // CSolver[iw].factorize( Csym );
  328. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  329. // }
  330. // }
  331. // template<>
  332. // void EMSchur3D< Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::NaturalOrdering<int> > > ::BuildCDirectSolver() {
  333. // CSolver = new Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::NaturalOrdering<int> > [Omegas.size()];
  334. // for (int iw=0; iw<Omegas.size(); ++iw) {
  335. // Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  336. // jsw_timer timer;
  337. // timer.begin();
  338. // /* Complex system */
  339. // std::cout << "CSymSimplicialLLT<NaturalOrdering> pattern analyzing C_" << iw << ",";
  340. // std::cout.flush();
  341. // CSolver[iw].analyzePattern( Csym );
  342. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  343. // /* factorize */
  344. // timer.begin();
  345. // std::cout << "CSymSimplicialLLT<NaturalOrdering> factorising C_" << iw << ", ";
  346. // std::cout.flush();
  347. // CSolver[iw].factorize( Csym );
  348. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  349. // }
  350. // }
  351. //
  352. // template<>
  353. // void EMSchur3D< Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > > ::BuildCDirectSolver() {
  354. // CSolver = new Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > [Omegas.size()];
  355. // for (int iw=0; iw<Omegas.size(); ++iw) {
  356. // //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  357. // jsw_timer timer;
  358. // timer.begin();
  359. // /* Complex system */
  360. // std::cout << "CSymSimplicialLLT<AMDOrdering> pattern analyzing C_" << iw << ",";
  361. // std::cout.flush();
  362. // CSolver[iw].analyzePattern( Cvec[iw] );
  363. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  364. // /* factorize */
  365. // timer.begin();
  366. // std::cout << "CSymSimplicialLLT<AMDOrdering> factorising C_" << iw << ", ";
  367. // std::cout.flush();
  368. // CSolver[iw].factorize( Cvec[iw] );
  369. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  370. // }
  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 -----