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

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