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

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811
  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. VectorXcr A = CSolver[iw].solve(Se);
  175. // // Solve Real system instead
  176. // The Real system is quasi-definite, though an LDLT decomposition exists, CHOLMOD doesn't find it.
  177. // An LU can be done on this, but compute performance is very similiar to the complex system, and diagonal pivoting
  178. // cannot be assumed to be best, hurting solve time.
  179. // /* EXPERIMENTAL */
  180. // VectorXr b2 = VectorXr::Zero(2*(unx+uny+unz));
  181. // b2.head(unx+uny+unz) = Se.real();
  182. // b2.tail(unx+uny+unz) = Se.imag();
  183. // VectorXr A2 = CReSolver[iw].solve(b2);
  184. // A.real() = A2.head( unx+uny+unz );
  185. // A.imag() = -A2.tail( unx+uny+unz ); // Due to decomp. negative!
  186. // /* END EXPERIMENTAL */
  187. VectorXcr ADiv = D*A; // ADiv == RHS == D C^I Se
  188. //VectorXcr Error = ((Cc.selfadjointView<Eigen::Lower>()*A).array() - Se.array());
  189. #if LOWER==1
  190. std::cout << "Using Eigen::Lower to calculate Error" << std::endl;
  191. VectorXcr Error = ((Cvec[iw].selfadjointView<Eigen::Lower>()*A).array() - Se.array());
  192. #elif UPPER==1
  193. std::cout << "Using Eigen::Upper to calculate Error" << std::endl;
  194. VectorXcr Error = ((Cvec[iw].selfadjointView<Eigen::Upper>()*A).array() - Se.array());
  195. #endif
  196. logio << "|| Div(A) || = " << ADiv.norm()
  197. << "\tInital solution error="<< Error.norm() // Iteritive info
  198. // << "\tSolver reported error="<< CSolver[iw].error() // Iteritive info
  199. << "\ttime " << timer.end() / 60. << " [m] "
  200. // << CSolver[iw].iterations() << " iterations"
  201. << std::endl;
  202. //VectorXcr ADivMAC = ADiv.array() * MAC.array().cast<Complex>();
  203. //logio << "|| Div(A) || on MAC grid " << ADivMAC.norm() << std::endl;
  204. /////////////////////
  205. // Solve for Phi
  206. logio << "Solving for Phi " << std::flush;
  207. timer.begin();
  208. tol = 1e-20;
  209. int success(2);
  210. success = implicitbicgstab(D, idx, ms, ADiv, Phi, CSolver[iw], max_it, tol, errorn, iter_done, logio);
  211. //Phi.array() *= MAC.array().cast<Complex>(); // remove phi from air regions
  212. /* Restart if necessary */
  213. int nrestart(1);
  214. // TODO send MAC to implicitbicgstab?
  215. while (success == 2 && nrestart < 18 && iter_done > 1) {
  216. success = implicitbicgstab(D, idx, ms, ADiv, Phi, CSolver[iw], max_it, tol, errorn, iter_done, logio);
  217. //Phi.array() *= MAC.array().cast<Complex>(); // remove phi from air regions
  218. nrestart += 1;
  219. }
  220. logio << "Implicit BiCGStab solution in " << iter_done << " iterations."
  221. << " with error " << std::setprecision(8) << std::scientific << errorn << std::endl;
  222. logio << "time "<< timer.end()/60. << " [m]" << std::endl;
  223. E = ms.array()*(D.transpose()*Phi).array(); // Temp, field due to charge
  224. /////////////////////////////////////
  225. // Compute A
  226. /////////////////////////////////////
  227. logio << "Solving for A using phi" << std::endl;
  228. std::cout << "Solving for A" << std::endl;
  229. max_it = nx*ny*nz;
  230. tol = 5e-16;
  231. errorn = 0;
  232. iter_done = 0;
  233. timer.begin();
  234. A = CSolver[iw].solve( (Se-E).eval() ); // UmfPack requires eval?
  235. VectorXcr ADiv2 = D*A;
  236. //logio << "|| Div(A) || = " << ADiv2.norm() ;
  237. //" in " << iter_done << " iterations"
  238. //<< " with error " << errorn << "\t";
  239. // Report error of solutions
  240. //Error = ((Cc.selfadjointView<Eigen::Lower>()*A).array() + E.array() - Se.array());
  241. #if LOWER==1
  242. std::cout << "Using Eigen::Lower to calculate Error" << std::endl;
  243. Error = ((Cvec[iw].selfadjointView<Eigen::Lower>()*A).array() + E.array() - Se.array());
  244. #elif UPPER==1
  245. std::cout << "Using Eigen::Upper to calculate Error" << std::endl;
  246. Error = ((Cvec[iw].selfadjointView<Eigen::Upper>()*A).array() + E.array() - Se.array());
  247. #endif
  248. // << "\tSolver reported error="<< CSolver[iw].error() // Iteritive info
  249. // << "\ttime " << timer.end() / 60. << " [m] "
  250. // << CSolver[iw].iterations() << " iterations"
  251. logio << "|| Div(A) || = " << ADiv2.norm()
  252. << "\tSolution error="<< Error.norm() // Iteritive info
  253. // << "\tSolver reported error="<< CSolver[iw].error() // Iteritive info
  254. << "\ttime " << timer.end() / 60. << " [m] "
  255. // << CSolver[iw].iterations() << " iterations"
  256. << std::endl << std::endl;
  257. logio << std::fixed << std::setprecision(2) << "\ttime " << timer.end()/60. << "\ttotal time " << timer2.end()/60. << std::endl;
  258. logio.close();
  259. //////////////////////////////////////
  260. // Update Fields and report
  261. E.array() = Complex(0,-Omegas[iw])*A.array() - (D.transpose()*Phi).array(); // Secondary Field Only
  262. VectorXcr B = StaggeredGridCurl(A);
  263. WriteVTKResults( ResFile+ to_string(isource), A, Se, E0, E , Phi, ADiv, ADiv2, B);
  264. //dpoint->Delete();
  265. return ;
  266. } // ----- end of method EMSchur3D::SolveSource -----
  267. //--------------------------------------------------------------------------------------
  268. // Class: EMSchur3DBase
  269. // Method: BuildCDirectSolver
  270. //--------------------------------------------------------------------------------------
  271. template < class Solver >
  272. void EMSchur3D<Solver>::BuildCDirectSolver ( ) {
  273. CSolver = new Solver[Omegas.size()];
  274. for (int iw=0; iw<Omegas.size(); ++iw) {
  275. jsw_timer timer;
  276. timer.begin();
  277. /* Complex system */
  278. /*
  279. std::cout << "Generic solver pattern analyzing C_" << iw << ",";
  280. std::cout.flush();
  281. CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  282. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  283. // factorize
  284. timer.begin();
  285. std::cout << "Generic solver factorising C_" << iw << ", ";
  286. std::cout.flush();
  287. CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  288. */
  289. std::cerr << "No solver Specified!" << iw << ",";
  290. exit(EXIT_FAILURE);
  291. //CSolver[iw].compute( Cvec[iw].selfadjointView< Eigen::Lower>() );
  292. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  293. }
  294. }
  295. #ifdef HAVE_SUPERLU
  296. template<>
  297. void EMSchur3D< Eigen::SuperLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::BuildCDirectSolver() {
  298. CSolver = new Eigen::SuperLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > [Omegas.size()];
  299. for (int iw=0; iw<Omegas.size(); ++iw) {
  300. jsw_timer timer;
  301. timer.begin();
  302. /* SuperLU */
  303. // Recommended values for symmetric mode
  304. CSolver[iw].options().DiagPivotThresh = 0.01;
  305. CSolver[iw].options().SymmetricMode = YES;
  306. CSolver[iw].options().ColPerm = MMD_AT_PLUS_A;
  307. //CSolver[iw].options().Trans = NOTRANS;
  308. //CSolver[iw].options().ConditionNumber = NO;
  309. //std::cout << "SuperLU options:\n";
  310. //std::cout << "\tPivot Threshold: " << CSolver[iw].options().DiagPivotThresh << std::endl;
  311. //std::cout << "\tSymmetric mode: " << CSolver[iw].options().SymmetricMode << std::endl;
  312. //std::cout << "\tEquilibrate: " << CSolver[iw].options().Equil << std::endl;
  313. //std::cout << "\tCol Permutation: " << CSolver[iw].options().ColPerm << std::endl;
  314. //std::cout << "\tTrans: " << CSolver[iw].options().Trans << std::endl;
  315. //std::cout << "\tCondition Number: " << CSolver[iw].options().ConditionNumber << std::endl;
  316. /* Complex system */
  317. std::cout << "SuperLU pattern analyzing C_" << iw << ",";
  318. std::cout.flush();
  319. CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  320. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  321. // factorize
  322. timer.begin();
  323. std::cout << "SuperLU factorising C_" << iw << ", ";
  324. std::cout.flush();
  325. CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  326. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  327. }
  328. }
  329. #endif
  330. template<>
  331. void EMSchur3D< Eigen::SparseLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::COLAMDOrdering<int> > >::BuildCDirectSolver() {
  332. CSolver = new Eigen::SparseLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::COLAMDOrdering<int> > [Omegas.size()];
  333. for (int iw=0; iw<Omegas.size(); ++iw) {
  334. jsw_timer timer;
  335. timer.begin();
  336. CSolver[iw].isSymmetric(true);
  337. //CSolver[iw].setPivotThreshold(0.0); // OK for symmetric complex systems with real and imaginary positive definite parts.
  338. // // but our imaginary part is negative definite
  339. //http://www.ams.org/journals/mcom/1998-67-224/S0025-5718-98-00978-8/S0025-5718-98-00978-8.pdf
  340. /* Complex system */
  341. std::cout << "SparseLU pattern analyzing C_" << iw << ",";
  342. std::cout.flush();
  343. CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  344. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  345. // factorize
  346. timer.begin();
  347. std::cout << "SparseLU factorising C_" << iw << ", ";
  348. std::cout.flush();
  349. CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  350. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  351. }
  352. }
  353. #ifdef HAVE_PARDISO
  354. template<>
  355. void EMSchur3D< Eigen::PardisoLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::BuildCDirectSolver() {
  356. CSolver = new Eigen::PardisoLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > [Omegas.size()];
  357. for (int iw=0; iw<Omegas.size(); ++iw) {
  358. jsw_timer timer;
  359. timer.begin();
  360. //CSolver[iw].pardisoParameterArray()[2] = 28; // OMP_NUM_THREADS?
  361. /* Complex system */
  362. std::cout << "PardisoLU pattern analyzing C_" << iw << ",";
  363. std::cout.flush();
  364. CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  365. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  366. // factorize
  367. timer.begin();
  368. std::cout << "PardisoLU factorising C_" << iw << ", ";
  369. std::cout.flush();
  370. CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  371. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  372. }
  373. }
  374. template<>
  375. void EMSchur3D< Eigen::PardisoLDLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Symmetric > >::BuildCDirectSolver() {
  376. CSolver = new Eigen::PardisoLDLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Symmetric > [Omegas.size()];
  377. for (int iw=0; iw<Omegas.size(); ++iw) {
  378. jsw_timer timer;
  379. timer.begin();
  380. //CSolver[iw].pardisoParameterArray()[2] = 28; // OMP_NUM_THREADS?
  381. //Eigen::SparseMatrix<Complex> Csym = Cvec[iw].selfadjointView< Eigen::Lower >();
  382. /* Complex system */
  383. std::cout << "PardisoLDLT pattern analyzing C_" << iw << ",";
  384. std::cout.flush();
  385. //CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  386. CSolver[iw].analyzePattern( Cvec[iw] );
  387. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  388. // factorize
  389. timer.begin();
  390. std::cout << "PardisoLDLT factorising C_" << iw << ", ";
  391. std::cout.flush();
  392. //CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  393. CSolver[iw].factorize( Cvec[iw] );
  394. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  395. }
  396. }
  397. #endif
  398. #ifdef HAVE_UMFPACK
  399. // Umfpack only seems to work when LOWER and UPPER are set to 1. Workarounds this have not been found.
  400. template<>
  401. void EMSchur3D< Eigen::UmfPackLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::BuildCDirectSolver() {
  402. CSolver = new Eigen::UmfPackLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > [Omegas.size()];
  403. for (int iw=0; iw<Omegas.size(); ++iw) {
  404. jsw_timer timer;
  405. timer.begin();
  406. /* Complex system */
  407. std::cout << "UmfPackLU pattern analyzing C_" << iw << ",";
  408. std::cout.flush();
  409. // Doesn't work, seg faults in solve
  410. //Eigen::SparseMatrix<Complex> Csym = Cvec[iw].selfadjointView< Eigen::Lower >();
  411. // Compiler errors get thrown with the view setting, good performance if LOWER and UPPER are set
  412. // CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() ); // Compiler error
  413. CSolver[iw].analyzePattern( Cvec[iw] ); // requires LOWER=1 UPPER=1, double memory
  414. // CSolver[iw].analyzePattern( Csym ); // seg faults
  415. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  416. // factorize
  417. timer.begin();
  418. std::cout << "UmfPackLU factorising C_" << iw << ", ";
  419. std::cout.flush();
  420. // CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  421. CSolver[iw].factorize( Cvec[iw] );
  422. // CSolver[iw].factorize( Csym );
  423. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  424. }
  425. }
  426. #endif
  427. // template<>
  428. // void EMSchur3D< Eigen::CholmodSupernodalLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
  429. // CSolver = new Eigen::CholmodSupernodalLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > [Omegas.size()];
  430. // for (int iw=0; iw<Omegas.size(); ++iw) {
  431. // Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  432. // jsw_timer timer;
  433. // timer.begin();
  434. // /* Complex system */
  435. // std::cout << "CholmodSupernodalLLT pattern analyzing C_" << iw << ",";
  436. // std::cout.flush();
  437. // CSolver[iw].analyzePattern( Csym );
  438. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  439. // /* factorize */
  440. // timer.begin();
  441. // std::cout << "CholmodSupernodalLLT factorising C_" << iw << ", ";
  442. // std::cout.flush();
  443. // CSolver[iw].factorize( Csym );
  444. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  445. // }
  446. // }
  447. // template<>
  448. // void EMSchur3D< Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::NaturalOrdering<int> > > ::BuildCDirectSolver() {
  449. // CSolver = new Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::NaturalOrdering<int> > [Omegas.size()];
  450. // for (int iw=0; iw<Omegas.size(); ++iw) {
  451. // Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  452. // jsw_timer timer;
  453. // timer.begin();
  454. // /* Complex system */
  455. // std::cout << "CSymSimplicialLLT<NaturalOrdering> pattern analyzing C_" << iw << ",";
  456. // std::cout.flush();
  457. // CSolver[iw].analyzePattern( Csym );
  458. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  459. // /* factorize */
  460. // timer.begin();
  461. // std::cout << "CSymSimplicialLLT<NaturalOrdering> factorising C_" << iw << ", ";
  462. // std::cout.flush();
  463. // CSolver[iw].factorize( Csym );
  464. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  465. // }
  466. // }
  467. //
  468. // template<>
  469. // void EMSchur3D< Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > > ::BuildCDirectSolver() {
  470. // CSolver = new Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > [Omegas.size()];
  471. // for (int iw=0; iw<Omegas.size(); ++iw) {
  472. // //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  473. // jsw_timer timer;
  474. // timer.begin();
  475. // /* Complex system */
  476. // std::cout << "CSymSimplicialLLT<AMDOrdering> pattern analyzing C_" << iw << ",";
  477. // std::cout.flush();
  478. // CSolver[iw].analyzePattern( Cvec[iw] );
  479. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  480. // /* factorize */
  481. // timer.begin();
  482. // std::cout << "CSymSimplicialLLT<AMDOrdering> factorising C_" << iw << ", ";
  483. // std::cout.flush();
  484. // CSolver[iw].factorize( Cvec[iw] );
  485. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  486. // }
  487. // }
  488. //
  489. // template<>
  490. // void EMSchur3D< Eigen::CSymSimplicialLDLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > > ::BuildCDirectSolver() {
  491. // CSolver = new Eigen::CSymSimplicialLDLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > [Omegas.size()];
  492. // for (int iw=0; iw<Omegas.size(); ++iw) {
  493. // Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  494. // jsw_timer timer;
  495. // timer.begin();
  496. // /* Complex system */
  497. // std::cout << "CSymSimplicialLDLT<AMDOrdering> pattern analyzing C_" << iw << ",";
  498. // std::cout.flush();
  499. // CSolver[iw].analyzePattern( Csym );
  500. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  501. // /* factorize */
  502. // timer.begin();
  503. // std::cout << "CSymSimplicialLDLT<AMDOrdering> factorising C_" << iw << ", ";
  504. // std::cout.flush();
  505. // CSolver[iw].factorize( Csym );
  506. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  507. // }
  508. // }
  509. template<>
  510. void EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::IncompleteLUT<Complex> > > ::BuildCDirectSolver() {
  511. CSolver = new Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::IncompleteLUT<Complex> > [Omegas.size()];
  512. for (int iw=0; iw<Omegas.size(); ++iw) {
  513. Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  514. CSolver[iw].preconditioner().setDroptol(1e-6); //1e-5); // 1e-12
  515. CSolver[iw].preconditioner().setFillfactor(5e1); // 1e2
  516. jsw_timer timer;
  517. timer.begin();
  518. /* Complex system */
  519. std::cout << "BiCGSTAB(ILU) pattern analyzing C_" << iw << ",";
  520. std::cout.flush();
  521. CSolver[iw].analyzePattern( Csym );
  522. //CSolver[iw].analyzePattern( Cvec[iw]);
  523. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  524. /* factorize */
  525. timer.begin();
  526. std::cout << "BiCGSTAB(ILU) factorising C_" << iw << ", ";
  527. std::cout.flush();
  528. CSolver[iw].factorize( Csym );
  529. //CSolver[iw].factorize( Cvec[iw] );
  530. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  531. }
  532. }
  533. template<>
  534. void EMSchur3D< Eigen::BiCGSTAB< Eigen::SparseMatrix<Complex, Eigen::ColMajor> > > ::BuildCDirectSolver() {
  535. CSolver = new Eigen::BiCGSTAB< Eigen::SparseMatrix<Complex, Eigen::ColMajor> > [Omegas.size()];
  536. for (int iw=0; iw<Omegas.size(); ++iw) {
  537. Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  538. jsw_timer timer;
  539. timer.begin();
  540. /* Complex system */
  541. std::cout << "BiCGSTAB pattern analyzing C_" << iw << ",";
  542. std::cout.flush();
  543. CSolver[iw].analyzePattern( Csym );
  544. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  545. // factorize
  546. timer.begin();
  547. std::cout << "BiCGSTAB factorising C_" << iw << ", ";
  548. std::cout.flush();
  549. CSolver[iw].factorize( Csym );
  550. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  551. }
  552. }
  553. template<>
  554. void EMSchur3D< Eigen::LeastSquaresConjugateGradient< Eigen::SparseMatrix<Complex, Eigen::ColMajor> > > ::BuildCDirectSolver() {
  555. CSolver = new Eigen::LeastSquaresConjugateGradient< Eigen::SparseMatrix<Complex, Eigen::ColMajor> > [Omegas.size()];
  556. for (int iw=0; iw<Omegas.size(); ++iw) {
  557. Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  558. jsw_timer timer;
  559. timer.begin();
  560. /* Complex system */
  561. std::cout << "LeastSquaresConjugateGradient pattern analyzing C_" << iw << ",";
  562. std::cout.flush();
  563. CSolver[iw].analyzePattern( Csym );
  564. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  565. // factorize
  566. timer.begin();
  567. std::cout << "LeastSquaresConjugateGradient factorising C_" << iw << ", ";
  568. std::cout.flush();
  569. CSolver[iw].factorize( Csym );
  570. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  571. }
  572. }
  573. template<>
  574. void EMSchur3D< Eigen::ConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
  575. CSolver = new Eigen::ConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > [Omegas.size()];
  576. for (int iw=0; iw<Omegas.size(); ++iw) {
  577. //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  578. jsw_timer timer;
  579. timer.begin();
  580. /* Complex system */
  581. std::cout << "ConjugateGradient pattern analyzing C_" << iw << ",";
  582. std::cout.flush();
  583. CSolver[iw].analyzePattern( Cvec[iw] );
  584. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  585. // factorize
  586. timer.begin();
  587. std::cout << "ConjugateGradient factorising C_" << iw << ", ";
  588. std::cout.flush();
  589. CSolver[iw].factorize( Cvec[iw] );
  590. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  591. }
  592. }
  593. template<>
  594. void EMSchur3D< Eigen::ConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::IncompleteCholesky<Complex> > > ::BuildCDirectSolver() {
  595. CSolver = new Eigen::ConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::IncompleteCholesky<Complex> > [Omegas.size()];
  596. for (int iw=0; iw<Omegas.size(); ++iw) {
  597. //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  598. jsw_timer timer;
  599. timer.begin();
  600. /* Complex system */
  601. std::cout << "ConjugateGradient<IncompleteCholesky> pattern analyzing C_" << iw << ",";
  602. std::cout.flush();
  603. CSolver[iw].analyzePattern( Cvec[iw] );
  604. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  605. // factorize
  606. timer.begin();
  607. std::cout << "ConjugateGradient<IncompleteCholesky> factorising C_" << iw << ", ";
  608. std::cout.flush();
  609. CSolver[iw].factorize( Cvec[iw] );
  610. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  611. }
  612. }
  613. // template<>
  614. // void EMSchur3D< Eigen::PastixLLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
  615. // CSolver = new Eigen::PastixLLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > [Omegas.size()];
  616. // //MPI_Init(NULL, NULL);
  617. // for (int iw=0; iw<Omegas.size(); ++iw) {
  618. // //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  619. // jsw_timer timer;
  620. // timer.begin();
  621. // /* Complex system */
  622. // std::cout << "PaStiX LLT pattern analyzing C_" << iw << ",";
  623. // std::cout.flush();
  624. // CSolver[iw].analyzePattern( Cvec[iw] );
  625. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  626. // // factorize
  627. // timer.begin();
  628. // std::cout << "PaStiX LLT factorising C_" << iw << ", ";
  629. // std::cout.flush();
  630. // CSolver[iw].factorize( Cvec[iw] );
  631. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  632. // }
  633. // }
  634. //
  635. // template<>
  636. // void EMSchur3D< Eigen::PastixLDLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
  637. // CSolver = new Eigen::PastixLDLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > [Omegas.size()];
  638. // //MPI_Init(NULL, NULL);
  639. // for (int iw=0; iw<Omegas.size(); ++iw) {
  640. // //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
  641. // jsw_timer timer;
  642. // timer.begin();
  643. // /* Complex system */
  644. // std::cout << "PaStiX LDLT pattern analyzing C_" << iw << ",";
  645. // std::cout.flush();
  646. // CSolver[iw].analyzePattern( Cvec[iw] );
  647. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  648. // // factorize
  649. // timer.begin();
  650. // std::cout << "PaStiX LDLT factorising C_" << iw << ", ";
  651. // std::cout.flush();
  652. // CSolver[iw].factorize( Cvec[iw] );
  653. // std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  654. // std::cout << "INFO " << CSolver[iw].info( ) << std::endl;
  655. // }
  656. // }
  657. //
  658. #ifdef HAVE_PASTIX
  659. template<>
  660. void EMSchur3D< Eigen::PastixLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, true > > ::BuildCDirectSolver() {
  661. CSolver = new Eigen::PastixLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, true > [Omegas.size()];
  662. for (int iw=0; iw<Omegas.size(); ++iw) {
  663. jsw_timer timer;
  664. timer.begin();
  665. /* Complex system */
  666. std::cout << "PastixLU pattern analyzing C_" << iw << ",";
  667. std::cout.flush();
  668. CSolver[iw].analyzePattern( Cvec[iw].selfadjointView< Eigen::Lower>() );
  669. //CSolver[iw].analyzePattern( Cvec[iw] );
  670. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  671. // factorize
  672. timer.begin();
  673. std::cout << "PastixLU factorising C_" << iw << ", ";
  674. std::cout.flush();
  675. CSolver[iw].factorize( Cvec[iw].selfadjointView< Eigen::Lower>() );
  676. //CSolver[iw].factorize( Cvec[iw] );
  677. std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
  678. }
  679. }
  680. #endif
  681. } // ----- end of Lemma name -----
  682. #endif // ----- #ifndef EMSCHUR3D_INC -----