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

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