3D EM based on Schur decomposition
Вы не можете выбрать более 25 тем Темы должны начинаться с буквы или цифры, могут содержать дефисы(-) и должны содержать не более 35 символов.

EMSchur3D.h 28KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655
  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] " << CSolver[iw].iterations() << " iterations" << std::endl;
  185. //VectorXcr ADivMAC = ADiv.array() * MAC.array().cast<Complex>();
  186. //logio << "|| Div(A) || on MAC grid " << ADivMAC.norm() << std::endl;
  187. /////////////////////
  188. // Solve for Phi
  189. logio << "Solving for Phi " << std::flush;
  190. timer.begin();
  191. tol = 1e-18;
  192. int success(2);
  193. success = implicitbicgstab(D, idx, ms, ADiv, Phi, CSolver[iw], max_it, tol, errorn, iter_done, logio);
  194. //Phi.array() *= MAC.array().cast<Complex>(); // remove phi from air regions
  195. /* Restart if necessary */
  196. /*
  197. int nrestart(1);
  198. // TODO send MAC to implicitbicgstab?
  199. while (success == 2 && nrestart < 18 && iter_done > 1) {
  200. success = implicitbicgstab(D, idx, ms, ADiv, Phi, CSolver[iw], max_it, tol, errorn, iter_done, logio);
  201. //Phi.array() *= MAC.array().cast<Complex>(); // remove phi from air regions
  202. nrestart += 1;
  203. }
  204. */
  205. logio << "Implicit BiCGStab solution in " << iter_done << " iterations."
  206. << " with error " << std::setprecision(8) << std::scientific << errorn << std::endl;
  207. logio << "time "<< timer.end()/60. << " [m]" << std::endl;
  208. E = ms.array()*(D.transpose()*Phi).array(); // Temp, field due to charge
  209. /////////////////////////////////////
  210. // Compute A
  211. /////////////////////////////////////
  212. logio << "Solving for A using phi" << std::endl;
  213. std::cout << "Solving for A" << std::endl;
  214. max_it = nx*ny*nz;
  215. tol = 5e-16;
  216. errorn = 0;
  217. iter_done = 0;
  218. timer.begin();
  219. A = CSolver[iw].solve( (Se-E).eval() ); // UmfPack requires eval?
  220. VectorXcr ADiv2 = D*A;
  221. logio << "|| Div(A) || = " << ADiv2.norm() ;
  222. //" in " << iter_done << " iterations"
  223. //<< " with error " << errorn << "\t";
  224. // Report error of solutions
  225. Error = ((Cc.selfadjointView<Eigen::Lower>()*A).array() + E.array() - Se.array());
  226. logio << "\tsolution error " << Error.norm()
  227. << std::fixed << std::setprecision(2) << "\ttime " << timer.end()/60. << "\ttotal time " << timer2.end()/60. << std::endl;
  228. logio.close();
  229. //////////////////////////////////////
  230. // Update Fields and report
  231. E.array() = Complex(0,-Omegas[iw])*A.array() - (D.transpose()*Phi).array(); // Secondary Field Only
  232. VectorXcr B = StaggeredGridCurl(A);
  233. WriteVTKResults( ResFile+ to_string(isource), A, Se, E0, E , Phi, ADiv, ADiv2, B);
  234. //dpoint->Delete();
  235. return ;
  236. } // ----- end of method EMSchur3D::SolveSource -----
  237. //--------------------------------------------------------------------------------------
  238. // Class: EMSchur3DBase
  239. // Method: BuildCDirectSolver
  240. //--------------------------------------------------------------------------------------
  241. template < class Solver >
  242. void EMSchur3D<Solver>::BuildCDirectSolver ( ) {
  243. CSolver = new Solver[Omegas.size()];
  244. for (int iw=0; iw<Omegas.size(); ++iw) {
  245. jsw_timer timer;
  246. timer.begin();
  247. /* Complex system */
  248. /*
  249. std::cout << "Generic solver pattern analyzing C_" << iw << ",";
  250. std::cout.flush();
  251. CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  252. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  253. // factorize
  254. timer.begin();
  255. std::cout << "Generic solver factorising C_" << iw << ", ";
  256. std::cout.flush();
  257. CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  258. */
  259. std::cerr << "No solver Specified!" << iw << ",";
  260. exit(EXIT_FAILURE);
  261. //CSolver[iw].compute( Cvec[iw].selfadjointView< Eigen::Lower>() );
  262. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  263. }
  264. }
  265. #ifdef HAVE_SUPERLUMT
  266. template<>
  267. void EMSchur3D< Eigen::SuperLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor> > >::BuildCDirectSolver() {
  268. CSolver = new Eigen::SuperLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor> > [Omegas.size()];
  269. for (int iw=0; iw<Omegas.size(); ++iw) {
  270. jsw_timer timer;
  271. timer.begin();
  272. /* SuperLU */
  273. //CSolver[iw].options().DiagPivotThresh = 0.01;
  274. //CSolver[iw].options().SymmetricMode = YES;
  275. //CSolver[iw].options().ColPerm = MMD_AT_PLUS_A;
  276. //CSolver[iw].options().Trans = NOTRANS;
  277. //CSolver[iw].options().ConditionNumber = NO;
  278. //std::cout << "SuperLU options:\n";
  279. //std::cout << "\tPivot Threshold: " << CSolver[iw].options().DiagPivotThresh << std::endl;
  280. //std::cout << "\tSymmetric mode: " << CSolver[iw].options().SymmetricMode << std::endl;
  281. //std::cout << "\tEquilibrate: " << CSolver[iw].options().Equil << std::endl;
  282. //std::cout << "\tCol Permutation: " << CSolver[iw].options().ColPerm << std::endl;
  283. //std::cout << "\tTrans: " << CSolver[iw].options().Trans << std::endl;
  284. //std::cout << "\tCondition Number: " << CSolver[iw].options().ConditionNumber << std::endl;
  285. /* Complex system */
  286. std::cout << "SuperLU_MT pattern analyzing C_" << iw << ",";
  287. std::cout.flush();
  288. CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  289. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  290. // factorize
  291. timer.begin();
  292. std::cout << "SuperLU_MT factorising C_" << iw << ", ";
  293. std::cout.flush();
  294. CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  295. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  296. }
  297. }
  298. #endif
  299. template<>
  300. void EMSchur3D< Eigen::SparseLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::COLAMDOrdering<int> > >::BuildCDirectSolver() {
  301. CSolver = new Eigen::SparseLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::COLAMDOrdering<int> > [Omegas.size()];
  302. for (int iw=0; iw<Omegas.size(); ++iw) {
  303. jsw_timer timer;
  304. timer.begin();
  305. CSolver[iw].isSymmetric(true);
  306. CSolver[iw].setPivotThreshold(0.0);
  307. /* Complex system */
  308. std::cout << "SparseLU pattern analyzing C_" << iw << ",";
  309. std::cout.flush();
  310. CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  311. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  312. // factorize
  313. timer.begin();
  314. std::cout << "SparseLU factorising C_" << iw << ", ";
  315. std::cout.flush();
  316. CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  317. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  318. }
  319. }
  320. // template<>
  321. // void EMSchur3D< Eigen::CholmodSupernodalLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
  322. // CSolver = new Eigen::CholmodSupernodalLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > [Omegas.size()];
  323. // for (int iw=0; iw<Omegas.size(); ++iw) {
  324. // Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  325. // jsw_timer timer;
  326. // timer.begin();
  327. // /* Complex system */
  328. // std::cout << "CholmodSupernodalLLT pattern analyzing C_" << iw << ",";
  329. // std::cout.flush();
  330. // CSolver[iw].analyzePattern( Csym );
  331. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  332. // /* factorize */
  333. // timer.begin();
  334. // std::cout << "CholmodSupernodalLLT factorising C_" << iw << ", ";
  335. // std::cout.flush();
  336. // CSolver[iw].factorize( Csym );
  337. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  338. // }
  339. // }
  340. // template<>
  341. // void EMSchur3D< Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::NaturalOrdering<int> > > ::BuildCDirectSolver() {
  342. // CSolver = new Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::NaturalOrdering<int> > [Omegas.size()];
  343. // for (int iw=0; iw<Omegas.size(); ++iw) {
  344. // Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  345. // jsw_timer timer;
  346. // timer.begin();
  347. // /* Complex system */
  348. // std::cout << "CSymSimplicialLLT<NaturalOrdering> pattern analyzing C_" << iw << ",";
  349. // std::cout.flush();
  350. // CSolver[iw].analyzePattern( Csym );
  351. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  352. // /* factorize */
  353. // timer.begin();
  354. // std::cout << "CSymSimplicialLLT<NaturalOrdering> factorising C_" << iw << ", ";
  355. // std::cout.flush();
  356. // CSolver[iw].factorize( Csym );
  357. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  358. // }
  359. // }
  360. //
  361. // template<>
  362. // void EMSchur3D< Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > > ::BuildCDirectSolver() {
  363. // CSolver = new Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > [Omegas.size()];
  364. // for (int iw=0; iw<Omegas.size(); ++iw) {
  365. // //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  366. // jsw_timer timer;
  367. // timer.begin();
  368. // /* Complex system */
  369. // std::cout << "CSymSimplicialLLT<AMDOrdering> pattern analyzing C_" << iw << ",";
  370. // std::cout.flush();
  371. // CSolver[iw].analyzePattern( Cvec[iw] );
  372. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  373. // /* factorize */
  374. // timer.begin();
  375. // std::cout << "CSymSimplicialLLT<AMDOrdering> factorising C_" << iw << ", ";
  376. // std::cout.flush();
  377. // CSolver[iw].factorize( Cvec[iw] );
  378. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  379. // }
  380. // }
  381. //
  382. // template<>
  383. // void EMSchur3D< Eigen::CSymSimplicialLDLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > > ::BuildCDirectSolver() {
  384. // CSolver = new Eigen::CSymSimplicialLDLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > [Omegas.size()];
  385. // for (int iw=0; iw<Omegas.size(); ++iw) {
  386. // Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  387. // jsw_timer timer;
  388. // timer.begin();
  389. // /* Complex system */
  390. // std::cout << "CSymSimplicialLDLT<AMDOrdering> pattern analyzing C_" << iw << ",";
  391. // std::cout.flush();
  392. // CSolver[iw].analyzePattern( Csym );
  393. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  394. // /* factorize */
  395. // timer.begin();
  396. // std::cout << "CSymSimplicialLDLT<AMDOrdering> factorising C_" << iw << ", ";
  397. // std::cout.flush();
  398. // CSolver[iw].factorize( Csym );
  399. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  400. // }
  401. // }
  402. template<>
  403. void EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::IncompleteLUT<Complex> > > ::BuildCDirectSolver() {
  404. CSolver = new Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::IncompleteLUT<Complex> > [Omegas.size()];
  405. for (int iw=0; iw<Omegas.size(); ++iw) {
  406. Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  407. CSolver[iw].preconditioner().setDroptol(1e-12); // 1e-12
  408. CSolver[iw].preconditioner().setFillfactor(1e1); // 1e2
  409. jsw_timer timer;
  410. timer.begin();
  411. /* Complex system */
  412. std::cout << "BiCGSTAB(ILU) pattern analyzing C_" << iw << ",";
  413. std::cout.flush();
  414. CSolver[iw].analyzePattern( Csym );
  415. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  416. /* factorize */
  417. timer.begin();
  418. std::cout << "BiCGSTAB(ILU) factorising C_" << iw << ", ";
  419. std::cout.flush();
  420. CSolver[iw].factorize( Csym );
  421. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  422. }
  423. }
  424. template<>
  425. void EMSchur3D< Eigen::BiCGSTAB< Eigen::SparseMatrix<Complex, Eigen::RowMajor> > > ::BuildCDirectSolver() {
  426. CSolver = new Eigen::BiCGSTAB< Eigen::SparseMatrix<Complex, Eigen::RowMajor> > [Omegas.size()];
  427. for (int iw=0; iw<Omegas.size(); ++iw) {
  428. Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  429. jsw_timer timer;
  430. timer.begin();
  431. /* Complex system */
  432. std::cout << "BiCGSTAB pattern analyzing C_" << iw << ",";
  433. std::cout.flush();
  434. CSolver[iw].analyzePattern( Csym );
  435. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  436. // factorize
  437. timer.begin();
  438. std::cout << "BiCGSTAB factorising C_" << iw << ", ";
  439. std::cout.flush();
  440. CSolver[iw].factorize( Csym );
  441. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  442. }
  443. }
  444. template<>
  445. void EMSchur3D< Eigen::LeastSquaresConjugateGradient< Eigen::SparseMatrix<Complex, Eigen::RowMajor> > > ::BuildCDirectSolver() {
  446. CSolver = new Eigen::LeastSquaresConjugateGradient< Eigen::SparseMatrix<Complex, Eigen::RowMajor> > [Omegas.size()];
  447. for (int iw=0; iw<Omegas.size(); ++iw) {
  448. Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  449. jsw_timer timer;
  450. timer.begin();
  451. /* Complex system */
  452. std::cout << "LeastSquaresConjugateGradient pattern analyzing C_" << iw << ",";
  453. std::cout.flush();
  454. CSolver[iw].analyzePattern( Csym );
  455. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  456. // factorize
  457. timer.begin();
  458. std::cout << "LeastSquaresConjugateGradient factorising C_" << iw << ", ";
  459. std::cout.flush();
  460. CSolver[iw].factorize( Csym );
  461. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  462. }
  463. }
  464. template<>
  465. void EMSchur3D< Eigen::ConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
  466. CSolver = new Eigen::ConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > [Omegas.size()];
  467. for (int iw=0; iw<Omegas.size(); ++iw) {
  468. //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  469. jsw_timer timer;
  470. timer.begin();
  471. /* Complex system */
  472. std::cout << "ConjugateGradient pattern analyzing C_" << iw << ",";
  473. std::cout.flush();
  474. CSolver[iw].analyzePattern( Cvec[iw] );
  475. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  476. // factorize
  477. timer.begin();
  478. std::cout << "ConjugateGradient factorising C_" << iw << ", ";
  479. std::cout.flush();
  480. CSolver[iw].factorize( Cvec[iw] );
  481. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  482. }
  483. }
  484. // template<>
  485. // void EMSchur3D< Eigen::PastixLLT<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
  486. // CSolver = new Eigen::PastixLLT<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > [Omegas.size()];
  487. // //MPI_Init(NULL, NULL);
  488. // for (int iw=0; iw<Omegas.size(); ++iw) {
  489. // //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  490. // jsw_timer timer;
  491. // timer.begin();
  492. // /* Complex system */
  493. // std::cout << "PaStiX LLT pattern analyzing C_" << iw << ",";
  494. // std::cout.flush();
  495. // CSolver[iw].analyzePattern( Cvec[iw] );
  496. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  497. // // factorize
  498. // timer.begin();
  499. // std::cout << "PaStiX LLT factorising C_" << iw << ", ";
  500. // std::cout.flush();
  501. // CSolver[iw].factorize( Cvec[iw] );
  502. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  503. // }
  504. // }
  505. //
  506. // template<>
  507. // void EMSchur3D< Eigen::PastixLDLT<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
  508. // CSolver = new Eigen::PastixLDLT<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > [Omegas.size()];
  509. // //MPI_Init(NULL, NULL);
  510. // for (int iw=0; iw<Omegas.size(); ++iw) {
  511. // //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  512. // jsw_timer timer;
  513. // timer.begin();
  514. // /* Complex system */
  515. // std::cout << "PaStiX LDLT pattern analyzing C_" << iw << ",";
  516. // std::cout.flush();
  517. // CSolver[iw].analyzePattern( Cvec[iw] );
  518. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  519. // // factorize
  520. // timer.begin();
  521. // std::cout << "PaStiX LDLT factorising C_" << iw << ", ";
  522. // std::cout.flush();
  523. // CSolver[iw].factorize( Cvec[iw] );
  524. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  525. // std::cout << "INFO " << CSolver[iw].info( ) << std::endl;
  526. // }
  527. // }
  528. //
  529. // template<>
  530. // void EMSchur3D< Eigen::PastixLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, true > > ::BuildCDirectSolver() {
  531. // CSolver = new Eigen::PastixLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, true > [Omegas.size()];
  532. // //MPI_Init(NULL, NULL);
  533. // for (int iw=0; iw<Omegas.size(); ++iw) {
  534. // Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  535. // jsw_timer timer;
  536. // timer.begin();
  537. // /* Complex system */
  538. // std::cout << "PaStiX LU pattern analyzing C_" << iw << ",";
  539. // std::cout.flush();
  540. // CSolver[iw].compute( Csym );
  541. // std::cout << "PaStiX LU Done C_" << iw << std::endl;;
  542. // // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  543. // // // factorize
  544. // // timer.begin();
  545. // // std::cout << "PaStiX LU factorising C_" << iw << ", ";
  546. // // std::cout.flush();
  547. // // CSolver[iw].factorize( Csym );
  548. // // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  549. // }
  550. // }
  551. } // ----- end of Lemma name -----
  552. #endif // ----- #ifndef EMSCHUR3D_INC -----