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

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680
  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"; // End of doc
  33. return stream;
  34. }
  35. //friend std::ostream &operator<<(std::ostream &stream,
  36. // const EMSchur3D &ob);
  37. public:
  38. // ==================== LIFECYCLE =======================
  39. /**
  40. * @copybrief LemmaObject::New()
  41. * @copydetails LemmaObject::New()
  42. */
  43. static std::shared_ptr< EMSchur3D > NewSP() {
  44. return std::make_shared< EMSchur3D<Solver> >( ctor_key() );
  45. //return std::make_shared< EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::RowMajor> > > >( ctor_key() ) ;
  46. }
  47. /** Default protected constructor, use New */
  48. explicit EMSchur3D ( const ctor_key& key ) : EMSchur3DBase( key ), CSolver( nullptr ) {
  49. }
  50. /** Locked DeDerializing constructor, use factory DeSerialize method*/
  51. EMSchur3D (const YAML::Node& node, const ctor_key& key): EMSchur3DBase(node, key), CSolver( nullptr ) {
  52. }
  53. /** Default protected destructor, use Delete */
  54. virtual ~EMSchur3D () {
  55. // TODO delete arrays
  56. }
  57. /**
  58. * Uses YAML to serialize this object.
  59. * @return a YAML::Node
  60. */
  61. YAML::Node Serialize() const {
  62. YAML::Node node = EMSchur3DBase::Serialize();
  63. //node["NumberOfLayers"] = NumberOfLayers;
  64. node.SetTag( this->GetName() );
  65. return node;
  66. }
  67. /**
  68. * Constructs an object from a YAML::Node.
  69. */
  70. static EMSchur3D* DeSerialize(const YAML::Node& node);
  71. // ==================== OPERATORS =======================
  72. // ==================== OPERATIONS =======================
  73. /** Solves a single source problem. This method is thread safe.
  74. * @param[in] Source is the source term for generating primary fields
  75. * @param[in] isource is the source index
  76. */
  77. void SolveSource( std::shared_ptr<DipoleSource> Source , const int& isource);
  78. /** Builds the solver for the C matrix */
  79. void BuildCDirectSolver( );
  80. // ==================== ACCESS =======================
  81. virtual std::string GetName() const {
  82. return this->CName;
  83. }
  84. // ==================== INQUIRY =======================
  85. protected:
  86. // ==================== LIFECYCLE =======================
  87. private:
  88. /** Copy constructor */
  89. EMSchur3D( const EMSchur3D& ) = delete;
  90. // ==================== DATA MEMBERS =========================
  91. /** The templated solver for C */
  92. Solver* CSolver;
  93. Eigen::SparseMatrix<Complex> Csym;
  94. static constexpr auto CName = "EMSchur3D";
  95. }; // ----- end of class EMSchur3D -----
  96. ////////////////////////////////////////////////////////////////////////////////////////
  97. // Implimentation and Specialisations //
  98. ////////////////////////////////////////////////////////////////////////////////////////
  99. //--------------------------------------------------------------------------------------
  100. // Class: EMSchur3D
  101. // Method: SolveSource
  102. //--------------------------------------------------------------------------------------
  103. template < class Solver >
  104. void EMSchur3D<Solver>::SolveSource ( std::shared_ptr<DipoleSource> Source, const int& isource ) {
  105. std::cout << "In vanilla SolveSource" << std::endl;
  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. std::cout << "Done with primary field" << std::endl;
  129. // Allocate a ton of memory
  130. VectorXcr Phi = VectorXcr::Zero(uns);
  131. VectorXcr ms(unx+uny+unz); // mu sigma
  132. // Vector potential (A) Vector and phi
  133. VectorXcr Se = VectorXcr::Zero(unx+uny+unz);
  134. //VectorXcr A = VectorXcr::Zero(unx+uny+unz);
  135. VectorXcr E = VectorXcr::Zero(unx+uny+unz);
  136. VectorXcr E0 = VectorXcr::Zero(unx+uny+unz);
  137. // Lets get cracking
  138. std::cout << "Filling source terms" << std::endl;
  139. FillSourceTerms(ms, Se, E0, dpoint, Omegas[iw]);
  140. std::cout << "Done source terms" << std::endl;
  141. /////////////////////////////////////////////////
  142. // LOG File
  143. std::string logfile (ResFile);
  144. logfile += to_string(isource) + std::string(".log");
  145. ofstream logio(logfile.c_str());
  146. std::cout << "just logging, TODO fix source" << std::endl;
  147. // logio << *Source << std::endl;
  148. logio << *Grid << std::endl;
  149. logio << *LayModel << std::endl;
  150. std::cout << "dun logging" << std::endl;
  151. // solve for RHS
  152. int max_it(nx*ny*nz), iter_done(0);
  153. Real tol(3e-16), errorn(0);
  154. logio << "solving RHS for source " << isource << std::endl;
  155. // TODO, this is stupid, try and get rid of this copy!
  156. Eigen::SparseMatrix<Complex> Cc = Cvec[iw];
  157. jsw_timer timer;
  158. jsw_timer timer2;
  159. timer.begin();
  160. timer2.begin();
  161. /////////////////////////////////////////
  162. // Solve for RHS
  163. //CSolver[iw].setMaxIterations(750);
  164. VectorXcr A = CSolver[iw].solve(Se);
  165. // // Solve Real system instead
  166. // The Real system is quasi-definite, though an LDLT decomposition exists, CHOLMOD doesn't find it.
  167. // An LU can be done on this, but compute performance is very similiar to the complex system, and diagonal pivoting
  168. // cannot be assumed to be best, hurting solve time.
  169. // /* EXPERIMENTAL */
  170. // VectorXr b2 = VectorXr::Zero(2*(unx+uny+unz));
  171. // b2.head(unx+uny+unz) = Se.real();
  172. // b2.tail(unx+uny+unz) = Se.imag();
  173. // VectorXr A2 = CReSolver[iw].solve(b2);
  174. // A.real() = A2.head( unx+uny+unz );
  175. // A.imag() = -A2.tail( unx+uny+unz ); // Due to decomp. negative!
  176. // /* END EXPERIMENTAL */
  177. VectorXcr ADiv = D*A; // ADiv == RHS == D C^I Se
  178. VectorXcr Error = ((Cc.selfadjointView<Eigen::Lower>()*A).array() - Se.array());
  179. logio << "|| Div(A) || = " << ADiv.norm()
  180. // << " in " << iter_done << " iterations"
  181. //<< " with error " << errorn << "\t"
  182. << "\tInital solution error="<< Error.norm() // Iteritive info
  183. // << "\tSolver reported error="<< CSolver[iw].error() // Iteritive info
  184. << "\ttime " << timer.end() / 60. << " [m] "
  185. // << CSolver[iw].iterations() << " iterations"
  186. << std::endl;
  187. //VectorXcr ADivMAC = ADiv.array() * MAC.array().cast<Complex>();
  188. //logio << "|| Div(A) || on MAC grid " << ADivMAC.norm() << std::endl;
  189. /////////////////////
  190. // Solve for Phi
  191. logio << "Solving for Phi " << std::flush;
  192. timer.begin();
  193. tol = 1e-18;
  194. int success(2);
  195. success = implicitbicgstab(D, idx, ms, ADiv, Phi, CSolver[iw], max_it, tol, errorn, iter_done, logio);
  196. //Phi.array() *= MAC.array().cast<Complex>(); // remove phi from air regions
  197. /* Restart if necessary */
  198. /*
  199. int nrestart(1);
  200. // TODO send MAC to implicitbicgstab?
  201. while (success == 2 && nrestart < 18 && iter_done > 1) {
  202. success = implicitbicgstab(D, idx, ms, ADiv, Phi, CSolver[iw], max_it, tol, errorn, iter_done, logio);
  203. //Phi.array() *= MAC.array().cast<Complex>(); // remove phi from air regions
  204. nrestart += 1;
  205. }
  206. */
  207. logio << "Implicit BiCGStab solution in " << iter_done << " iterations."
  208. << " with error " << std::setprecision(8) << std::scientific << errorn << std::endl;
  209. logio << "time "<< timer.end()/60. << " [m]" << std::endl;
  210. E = ms.array()*(D.transpose()*Phi).array(); // Temp, field due to charge
  211. /////////////////////////////////////
  212. // Compute A
  213. /////////////////////////////////////
  214. logio << "Solving for A using phi" << std::endl;
  215. std::cout << "Solving for A" << std::endl;
  216. max_it = nx*ny*nz;
  217. tol = 5e-16;
  218. errorn = 0;
  219. iter_done = 0;
  220. timer.begin();
  221. A = CSolver[iw].solve( (Se-E).eval() ); // UmfPack requires eval?
  222. VectorXcr ADiv2 = D*A;
  223. logio << "|| Div(A) || = " << ADiv2.norm() ;
  224. //" in " << iter_done << " iterations"
  225. //<< " with error " << errorn << "\t";
  226. // Report error of solutions
  227. Error = ((Cc.selfadjointView<Eigen::Lower>()*A).array() + E.array() - Se.array());
  228. logio << "\tsolution error " << Error.norm()
  229. << std::fixed << std::setprecision(2) << "\ttime " << timer.end()/60. << "\ttotal time " << timer2.end()/60. << std::endl;
  230. logio.close();
  231. //////////////////////////////////////
  232. // Update Fields and report
  233. E.array() = Complex(0,-Omegas[iw])*A.array() - (D.transpose()*Phi).array(); // Secondary Field Only
  234. VectorXcr B = StaggeredGridCurl(A);
  235. WriteVTKResults( ResFile+ to_string(isource), A, Se, E0, E , Phi, ADiv, ADiv2, B);
  236. //dpoint->Delete();
  237. return ;
  238. } // ----- end of method EMSchur3D::SolveSource -----
  239. //--------------------------------------------------------------------------------------
  240. // Class: EMSchur3DBase
  241. // Method: BuildCDirectSolver
  242. //--------------------------------------------------------------------------------------
  243. template < class Solver >
  244. void EMSchur3D<Solver>::BuildCDirectSolver ( ) {
  245. CSolver = new Solver[Omegas.size()];
  246. for (int iw=0; iw<Omegas.size(); ++iw) {
  247. jsw_timer timer;
  248. timer.begin();
  249. /* Complex system */
  250. /*
  251. std::cout << "Generic solver pattern analyzing C_" << iw << ",";
  252. std::cout.flush();
  253. CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  254. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  255. // factorize
  256. timer.begin();
  257. std::cout << "Generic solver factorising C_" << iw << ", ";
  258. std::cout.flush();
  259. CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  260. */
  261. std::cerr << "No solver Specified!" << iw << ",";
  262. exit(EXIT_FAILURE);
  263. //CSolver[iw].compute( Cvec[iw].selfadjointView< Eigen::Lower>() );
  264. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  265. }
  266. }
  267. #ifdef HAVE_SUPERLUMT
  268. template<>
  269. void EMSchur3D< Eigen::SuperLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor> > >::BuildCDirectSolver() {
  270. CSolver = new Eigen::SuperLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor> > [Omegas.size()];
  271. for (int iw=0; iw<Omegas.size(); ++iw) {
  272. jsw_timer timer;
  273. timer.begin();
  274. /* SuperLU */
  275. //CSolver[iw].options().DiagPivotThresh = 0.01;
  276. //CSolver[iw].options().SymmetricMode = YES;
  277. //CSolver[iw].options().ColPerm = MMD_AT_PLUS_A;
  278. //CSolver[iw].options().Trans = NOTRANS;
  279. //CSolver[iw].options().ConditionNumber = NO;
  280. //std::cout << "SuperLU options:\n";
  281. //std::cout << "\tPivot Threshold: " << CSolver[iw].options().DiagPivotThresh << std::endl;
  282. //std::cout << "\tSymmetric mode: " << CSolver[iw].options().SymmetricMode << std::endl;
  283. //std::cout << "\tEquilibrate: " << CSolver[iw].options().Equil << std::endl;
  284. //std::cout << "\tCol Permutation: " << CSolver[iw].options().ColPerm << std::endl;
  285. //std::cout << "\tTrans: " << CSolver[iw].options().Trans << std::endl;
  286. //std::cout << "\tCondition Number: " << CSolver[iw].options().ConditionNumber << std::endl;
  287. /* Complex system */
  288. std::cout << "SuperLU_MT pattern analyzing C_" << iw << ",";
  289. std::cout.flush();
  290. CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  291. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  292. // factorize
  293. timer.begin();
  294. std::cout << "SuperLU_MT factorising C_" << iw << ", ";
  295. std::cout.flush();
  296. CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  297. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  298. }
  299. }
  300. #endif
  301. template<>
  302. void EMSchur3D< Eigen::SparseLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::COLAMDOrdering<int> > >::BuildCDirectSolver() {
  303. CSolver = new Eigen::SparseLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::COLAMDOrdering<int> > [Omegas.size()];
  304. for (int iw=0; iw<Omegas.size(); ++iw) {
  305. jsw_timer timer;
  306. timer.begin();
  307. CSolver[iw].isSymmetric(true);
  308. //CSolver[iw].setPivotThreshold(0.0); // OK for symmetric complex systems with real and imaginary positive definite parts.
  309. // // but our imaginary part is negative definite
  310. //http://www.ams.org/journals/mcom/1998-67-224/S0025-5718-98-00978-8/S0025-5718-98-00978-8.pdf
  311. /* Complex system */
  312. std::cout << "SparseLU pattern analyzing C_" << iw << ",";
  313. std::cout.flush();
  314. CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  315. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  316. // factorize
  317. timer.begin();
  318. std::cout << "SparseLU factorising C_" << iw << ", ";
  319. std::cout.flush();
  320. CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  321. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  322. }
  323. }
  324. // template<>
  325. // void EMSchur3D< Eigen::CholmodSupernodalLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
  326. // CSolver = new Eigen::CholmodSupernodalLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > [Omegas.size()];
  327. // for (int iw=0; iw<Omegas.size(); ++iw) {
  328. // Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  329. // jsw_timer timer;
  330. // timer.begin();
  331. // /* Complex system */
  332. // std::cout << "CholmodSupernodalLLT pattern analyzing C_" << iw << ",";
  333. // std::cout.flush();
  334. // CSolver[iw].analyzePattern( Csym );
  335. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  336. // /* factorize */
  337. // timer.begin();
  338. // std::cout << "CholmodSupernodalLLT factorising C_" << iw << ", ";
  339. // std::cout.flush();
  340. // CSolver[iw].factorize( Csym );
  341. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  342. // }
  343. // }
  344. // template<>
  345. // void EMSchur3D< Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::NaturalOrdering<int> > > ::BuildCDirectSolver() {
  346. // CSolver = new Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::NaturalOrdering<int> > [Omegas.size()];
  347. // for (int iw=0; iw<Omegas.size(); ++iw) {
  348. // Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  349. // jsw_timer timer;
  350. // timer.begin();
  351. // /* Complex system */
  352. // std::cout << "CSymSimplicialLLT<NaturalOrdering> pattern analyzing C_" << iw << ",";
  353. // std::cout.flush();
  354. // CSolver[iw].analyzePattern( Csym );
  355. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  356. // /* factorize */
  357. // timer.begin();
  358. // std::cout << "CSymSimplicialLLT<NaturalOrdering> factorising C_" << iw << ", ";
  359. // std::cout.flush();
  360. // CSolver[iw].factorize( Csym );
  361. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  362. // }
  363. // }
  364. //
  365. // template<>
  366. // void EMSchur3D< Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > > ::BuildCDirectSolver() {
  367. // CSolver = new Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > [Omegas.size()];
  368. // for (int iw=0; iw<Omegas.size(); ++iw) {
  369. // //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  370. // jsw_timer timer;
  371. // timer.begin();
  372. // /* Complex system */
  373. // std::cout << "CSymSimplicialLLT<AMDOrdering> pattern analyzing C_" << iw << ",";
  374. // std::cout.flush();
  375. // CSolver[iw].analyzePattern( Cvec[iw] );
  376. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  377. // /* factorize */
  378. // timer.begin();
  379. // std::cout << "CSymSimplicialLLT<AMDOrdering> factorising C_" << iw << ", ";
  380. // std::cout.flush();
  381. // CSolver[iw].factorize( Cvec[iw] );
  382. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  383. // }
  384. // }
  385. //
  386. // template<>
  387. // void EMSchur3D< Eigen::CSymSimplicialLDLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > > ::BuildCDirectSolver() {
  388. // CSolver = new Eigen::CSymSimplicialLDLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > [Omegas.size()];
  389. // for (int iw=0; iw<Omegas.size(); ++iw) {
  390. // Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  391. // jsw_timer timer;
  392. // timer.begin();
  393. // /* Complex system */
  394. // std::cout << "CSymSimplicialLDLT<AMDOrdering> pattern analyzing C_" << iw << ",";
  395. // std::cout.flush();
  396. // CSolver[iw].analyzePattern( Csym );
  397. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  398. // /* factorize */
  399. // timer.begin();
  400. // std::cout << "CSymSimplicialLDLT<AMDOrdering> factorising C_" << iw << ", ";
  401. // std::cout.flush();
  402. // CSolver[iw].factorize( Csym );
  403. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  404. // }
  405. // }
  406. template<>
  407. void EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::IncompleteLUT<Complex> > > ::BuildCDirectSolver() {
  408. CSolver = new Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::IncompleteLUT<Complex> > [Omegas.size()];
  409. for (int iw=0; iw<Omegas.size(); ++iw) {
  410. Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  411. CSolver[iw].preconditioner().setDroptol(1e-12); // 1e-12
  412. CSolver[iw].preconditioner().setFillfactor(1e1); // 1e2
  413. jsw_timer timer;
  414. timer.begin();
  415. /* Complex system */
  416. std::cout << "BiCGSTAB(ILU) pattern analyzing C_" << iw << ",";
  417. std::cout.flush();
  418. CSolver[iw].analyzePattern( Csym );
  419. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  420. /* factorize */
  421. timer.begin();
  422. std::cout << "BiCGSTAB(ILU) factorising C_" << iw << ", ";
  423. std::cout.flush();
  424. CSolver[iw].factorize( Csym );
  425. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  426. }
  427. }
  428. template<>
  429. void EMSchur3D< Eigen::BiCGSTAB< Eigen::SparseMatrix<Complex, Eigen::RowMajor> > > ::BuildCDirectSolver() {
  430. CSolver = new Eigen::BiCGSTAB< Eigen::SparseMatrix<Complex, Eigen::RowMajor> > [Omegas.size()];
  431. for (int iw=0; iw<Omegas.size(); ++iw) {
  432. Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  433. jsw_timer timer;
  434. timer.begin();
  435. /* Complex system */
  436. std::cout << "BiCGSTAB pattern analyzing C_" << iw << ",";
  437. std::cout.flush();
  438. CSolver[iw].analyzePattern( Csym );
  439. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  440. // factorize
  441. timer.begin();
  442. std::cout << "BiCGSTAB factorising C_" << iw << ", ";
  443. std::cout.flush();
  444. CSolver[iw].factorize( Csym );
  445. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  446. }
  447. }
  448. template<>
  449. void EMSchur3D< Eigen::LeastSquaresConjugateGradient< Eigen::SparseMatrix<Complex, Eigen::RowMajor> > > ::BuildCDirectSolver() {
  450. CSolver = new Eigen::LeastSquaresConjugateGradient< Eigen::SparseMatrix<Complex, Eigen::RowMajor> > [Omegas.size()];
  451. for (int iw=0; iw<Omegas.size(); ++iw) {
  452. Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  453. jsw_timer timer;
  454. timer.begin();
  455. /* Complex system */
  456. std::cout << "LeastSquaresConjugateGradient pattern analyzing C_" << iw << ",";
  457. std::cout.flush();
  458. CSolver[iw].analyzePattern( Csym );
  459. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  460. // factorize
  461. timer.begin();
  462. std::cout << "LeastSquaresConjugateGradient factorising C_" << iw << ", ";
  463. std::cout.flush();
  464. CSolver[iw].factorize( Csym );
  465. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  466. }
  467. }
  468. template<>
  469. void EMSchur3D< Eigen::ConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
  470. CSolver = new Eigen::ConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > [Omegas.size()];
  471. for (int iw=0; iw<Omegas.size(); ++iw) {
  472. //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  473. jsw_timer timer;
  474. timer.begin();
  475. /* Complex system */
  476. std::cout << "ConjugateGradient pattern analyzing C_" << iw << ",";
  477. std::cout.flush();
  478. CSolver[iw].analyzePattern( Cvec[iw] );
  479. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  480. // factorize
  481. timer.begin();
  482. std::cout << "ConjugateGradient factorising C_" << iw << ", ";
  483. std::cout.flush();
  484. CSolver[iw].factorize( Cvec[iw] );
  485. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  486. }
  487. }
  488. template<>
  489. void EMSchur3D< Eigen::ConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::IncompleteCholesky<Complex> > > ::BuildCDirectSolver() {
  490. CSolver = new Eigen::ConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::IncompleteCholesky<Complex> > [Omegas.size()];
  491. for (int iw=0; iw<Omegas.size(); ++iw) {
  492. //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  493. jsw_timer timer;
  494. timer.begin();
  495. /* Complex system */
  496. std::cout << "ConjugateGradient<IncompleteCholesky> pattern analyzing C_" << iw << ",";
  497. std::cout.flush();
  498. CSolver[iw].analyzePattern( Cvec[iw] );
  499. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  500. // factorize
  501. timer.begin();
  502. std::cout << "ConjugateGradient<IncompleteCholesky> factorising C_" << iw << ", ";
  503. std::cout.flush();
  504. CSolver[iw].factorize( Cvec[iw] );
  505. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  506. }
  507. }
  508. // template<>
  509. // void EMSchur3D< Eigen::PastixLLT<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
  510. // CSolver = new Eigen::PastixLLT<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > [Omegas.size()];
  511. // //MPI_Init(NULL, NULL);
  512. // for (int iw=0; iw<Omegas.size(); ++iw) {
  513. // //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  514. // jsw_timer timer;
  515. // timer.begin();
  516. // /* Complex system */
  517. // std::cout << "PaStiX LLT pattern analyzing C_" << iw << ",";
  518. // std::cout.flush();
  519. // CSolver[iw].analyzePattern( Cvec[iw] );
  520. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  521. // // factorize
  522. // timer.begin();
  523. // std::cout << "PaStiX LLT factorising C_" << iw << ", ";
  524. // std::cout.flush();
  525. // CSolver[iw].factorize( Cvec[iw] );
  526. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  527. // }
  528. // }
  529. //
  530. // template<>
  531. // void EMSchur3D< Eigen::PastixLDLT<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
  532. // CSolver = new Eigen::PastixLDLT<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > [Omegas.size()];
  533. // //MPI_Init(NULL, NULL);
  534. // for (int iw=0; iw<Omegas.size(); ++iw) {
  535. // //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  536. // jsw_timer timer;
  537. // timer.begin();
  538. // /* Complex system */
  539. // std::cout << "PaStiX LDLT pattern analyzing C_" << iw << ",";
  540. // std::cout.flush();
  541. // CSolver[iw].analyzePattern( Cvec[iw] );
  542. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  543. // // factorize
  544. // timer.begin();
  545. // std::cout << "PaStiX LDLT factorising C_" << iw << ", ";
  546. // std::cout.flush();
  547. // CSolver[iw].factorize( Cvec[iw] );
  548. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  549. // std::cout << "INFO " << CSolver[iw].info( ) << std::endl;
  550. // }
  551. // }
  552. //
  553. // template<>
  554. // void EMSchur3D< Eigen::PastixLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, true > > ::BuildCDirectSolver() {
  555. // CSolver = new Eigen::PastixLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, true > [Omegas.size()];
  556. // //MPI_Init(NULL, NULL);
  557. // for (int iw=0; iw<Omegas.size(); ++iw) {
  558. // Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  559. // jsw_timer timer;
  560. // timer.begin();
  561. // /* Complex system */
  562. // std::cout << "PaStiX LU pattern analyzing C_" << iw << ",";
  563. // std::cout.flush();
  564. // CSolver[iw].compute( Csym );
  565. // std::cout << "PaStiX LU Done C_" << iw << std::endl;;
  566. // // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  567. // // // factorize
  568. // // timer.begin();
  569. // // std::cout << "PaStiX LU factorising C_" << iw << ", ";
  570. // // std::cout.flush();
  571. // // CSolver[iw].factorize( Csym );
  572. // // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  573. // }
  574. // }
  575. } // ----- end of Lemma name -----
  576. #endif // ----- #ifndef EMSCHUR3D_INC -----