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

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