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

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