3D EM based on Schur decomposition
Vous ne pouvez pas sélectionner plus de 25 sujets Les noms de sujets doivent commencer par une lettre ou un nombre, peuvent contenir des tirets ('-') et peuvent comporter jusqu'à 35 caractères.

EMSchur3D.h 34KB

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