sp_auxlib_meat.hpp 48 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554
  1. // Copyright 2008-2016 Conrad Sanderson (http://conradsanderson.id.au)
  2. // Copyright 2008-2016 National ICT Australia (NICTA)
  3. //
  4. // Licensed under the Apache License, Version 2.0 (the "License");
  5. // you may not use this file except in compliance with the License.
  6. // You may obtain a copy of the License at
  7. // http://www.apache.org/licenses/LICENSE-2.0
  8. //
  9. // Unless required by applicable law or agreed to in writing, software
  10. // distributed under the License is distributed on an "AS IS" BASIS,
  11. // WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
  12. // See the License for the specific language governing permissions and
  13. // limitations under the License.
  14. // ------------------------------------------------------------------------
  15. //! \addtogroup sp_auxlib
  16. //! @{
  17. inline
  18. sp_auxlib::form_type
  19. sp_auxlib::interpret_form_str(const char* form_str)
  20. {
  21. arma_extra_debug_sigprint();
  22. // the order of the 3 if statements below is important
  23. if( form_str == NULL ) { return form_none; }
  24. if( form_str[0] == char(0) ) { return form_none; }
  25. if( form_str[1] == char(0) ) { return form_none; }
  26. const char c1 = form_str[0];
  27. const char c2 = form_str[1];
  28. if(c1 == 'l')
  29. {
  30. if(c2 == 'm') { return form_lm; }
  31. if(c2 == 'r') { return form_lr; }
  32. if(c2 == 'i') { return form_li; }
  33. if(c2 == 'a') { return form_la; }
  34. }
  35. else
  36. if(c1 == 's')
  37. {
  38. if(c2 == 'm') { return form_sm; }
  39. if(c2 == 'r') { return form_sr; }
  40. if(c2 == 'i') { return form_si; }
  41. if(c2 == 'a') { return form_sa; }
  42. }
  43. return form_none;
  44. }
  45. //! immediate eigendecomposition of symmetric real sparse object
  46. template<typename eT, typename T1>
  47. inline
  48. bool
  49. sp_auxlib::eigs_sym(Col<eT>& eigval, Mat<eT>& eigvec, const SpBase<eT, T1>& X, const uword n_eigvals, const char* form_str, const eT default_tol)
  50. {
  51. arma_extra_debug_sigprint();
  52. const unwrap_spmat<T1> U(X.get_ref());
  53. if((arma_config::debug) && (sp_auxlib::rudimentary_sym_check(U.M) == false))
  54. {
  55. if(is_cx<eT>::no ) { arma_debug_warn("eigs_sym(): given matrix is not symmetric"); }
  56. if(is_cx<eT>::yes) { arma_debug_warn("eigs_sym(): given matrix is not hermitian"); }
  57. }
  58. #if defined(ARMA_USE_NEWARP)
  59. {
  60. return sp_auxlib::eigs_sym_newarp(eigval, eigvec, U.M, n_eigvals, form_str, default_tol);
  61. }
  62. #elif defined(ARMA_USE_ARPACK)
  63. {
  64. return sp_auxlib::eigs_sym_arpack(eigval, eigvec, U.M, n_eigvals, form_str, default_tol);
  65. }
  66. #else
  67. {
  68. arma_ignore(eigval);
  69. arma_ignore(eigvec);
  70. arma_ignore(n_eigvals);
  71. arma_ignore(form_str);
  72. arma_ignore(default_tol);
  73. arma_stop_logic_error("eigs_sym(): use of NEWARP or ARPACK must be enabled");
  74. return false;
  75. }
  76. #endif
  77. }
  78. template<typename eT>
  79. inline
  80. bool
  81. sp_auxlib::eigs_sym_newarp(Col<eT>& eigval, Mat<eT>& eigvec, const SpMat<eT>& X, const uword n_eigvals, const char* form_str, const eT default_tol)
  82. {
  83. arma_extra_debug_sigprint();
  84. #if defined(ARMA_USE_NEWARP)
  85. {
  86. const form_type form_val = sp_auxlib::interpret_form_str(form_str);
  87. arma_debug_check( (form_val != form_lm) && (form_val != form_sm) && (form_val != form_la) && (form_val != form_sa), "eigs_sym(): unknown form specified" );
  88. const newarp::SparseGenMatProd<eT> op(X);
  89. arma_debug_check( (op.n_rows != op.n_cols), "eigs_sym(): given matrix must be square sized" );
  90. arma_debug_check( (n_eigvals >= op.n_rows), "eigs_sym(): n_eigvals must be less than the number of rows in the matrix" );
  91. // If the matrix is empty, the case is trivial.
  92. if( (op.n_cols == 0) || (n_eigvals == 0) ) // We already know n_cols == n_rows.
  93. {
  94. eigval.reset();
  95. eigvec.reset();
  96. return true;
  97. }
  98. uword n = op.n_rows;
  99. uword ncv = n_eigvals + 2 + 1;
  100. if(ncv < (2 * n_eigvals + 1)) { ncv = 2 * n_eigvals + 1; }
  101. if(ncv > n) { ncv = n; }
  102. eT tol = (std::max)(default_tol, std::numeric_limits<eT>::epsilon());
  103. // eigval.set_size(n_eigvals);
  104. // eigvec.set_size(n, n_eigvals);
  105. bool status = true;
  106. uword nconv = 0;
  107. try
  108. {
  109. if(form_val == form_lm)
  110. {
  111. newarp::SymEigsSolver< eT, newarp::EigsSelect::LARGEST_MAGN, newarp::SparseGenMatProd<eT> > eigs(op, n_eigvals, ncv);
  112. eigs.init();
  113. nconv = eigs.compute(1000, tol);
  114. eigval = eigs.eigenvalues();
  115. eigvec = eigs.eigenvectors();
  116. }
  117. else
  118. if(form_val == form_sm)
  119. {
  120. newarp::SymEigsSolver< eT, newarp::EigsSelect::SMALLEST_MAGN, newarp::SparseGenMatProd<eT> > eigs(op, n_eigvals, ncv);
  121. eigs.init();
  122. nconv = eigs.compute(1000, tol);
  123. eigval = eigs.eigenvalues();
  124. eigvec = eigs.eigenvectors();
  125. }
  126. else
  127. if(form_val == form_la)
  128. {
  129. newarp::SymEigsSolver< eT, newarp::EigsSelect::LARGEST_ALGE, newarp::SparseGenMatProd<eT> > eigs(op, n_eigvals, ncv);
  130. eigs.init();
  131. nconv = eigs.compute(1000, tol);
  132. eigval = eigs.eigenvalues();
  133. eigvec = eigs.eigenvectors();
  134. }
  135. else
  136. if(form_val == form_sa)
  137. {
  138. newarp::SymEigsSolver< eT, newarp::EigsSelect::SMALLEST_ALGE, newarp::SparseGenMatProd<eT> > eigs(op, n_eigvals, ncv);
  139. eigs.init();
  140. nconv = eigs.compute(1000, tol);
  141. eigval = eigs.eigenvalues();
  142. eigvec = eigs.eigenvectors();
  143. }
  144. }
  145. catch(const std::runtime_error&)
  146. {
  147. status = false;
  148. }
  149. if(status == true)
  150. {
  151. if(nconv == 0) { status = false; }
  152. }
  153. return status;
  154. }
  155. #else
  156. {
  157. arma_ignore(eigval);
  158. arma_ignore(eigvec);
  159. arma_ignore(X);
  160. arma_ignore(n_eigvals);
  161. arma_ignore(form_str);
  162. arma_ignore(default_tol);
  163. return false;
  164. }
  165. #endif
  166. }
  167. template<typename eT>
  168. inline
  169. bool
  170. sp_auxlib::eigs_sym_arpack(Col<eT>& eigval, Mat<eT>& eigvec, const SpMat<eT>& X, const uword n_eigvals, const char* form_str, const eT default_tol)
  171. {
  172. arma_extra_debug_sigprint();
  173. #if defined(ARMA_USE_ARPACK)
  174. {
  175. const form_type form_val = sp_auxlib::interpret_form_str(form_str);
  176. arma_debug_check( (form_val != form_lm) && (form_val != form_sm) && (form_val != form_la) && (form_val != form_sa), "eigs_sym(): unknown form specified" );
  177. char which_sm[3] = "SM";
  178. char which_lm[3] = "LM";
  179. char which_sa[3] = "SA";
  180. char which_la[3] = "LA";
  181. char* which;
  182. switch (form_val)
  183. {
  184. case form_sm: which = which_sm; break;
  185. case form_lm: which = which_lm; break;
  186. case form_sa: which = which_sa; break;
  187. case form_la: which = which_la; break;
  188. default: which = which_lm; break;
  189. }
  190. // Make sure it's square.
  191. arma_debug_check( (X.n_rows != X.n_cols), "eigs_sym(): given matrix must be square sized" );
  192. // Make sure we aren't asking for every eigenvalue.
  193. // The _saupd() functions allow asking for one more eigenvalue than the _naupd() functions.
  194. arma_debug_check( (n_eigvals >= X.n_rows), "eigs_sym(): n_eigvals must be less than the number of rows in the matrix" );
  195. // If the matrix is empty, the case is trivial.
  196. if( (X.n_cols == 0) || (n_eigvals == 0) ) // We already know n_cols == n_rows.
  197. {
  198. eigval.reset();
  199. eigvec.reset();
  200. return true;
  201. }
  202. // Set up variables that get used for neupd().
  203. blas_int n, ncv, ldv, lworkl, info;
  204. eT tol = default_tol;
  205. podarray<eT> resid, v, workd, workl;
  206. podarray<blas_int> iparam, ipntr;
  207. podarray<eT> rwork; // Not used in this case.
  208. run_aupd(n_eigvals, which, X, true /* sym, not gen */, n, tol, resid, ncv, v, ldv, iparam, ipntr, workd, workl, lworkl, rwork, info);
  209. if(info != 0)
  210. {
  211. return false;
  212. }
  213. // The process has converged, and now we need to recover the actual eigenvectors using seupd()
  214. blas_int rvec = 1; // .TRUE
  215. blas_int nev = n_eigvals;
  216. char howmny = 'A';
  217. char bmat = 'I'; // We are considering the standard eigenvalue problem.
  218. podarray<blas_int> select(ncv); // Logical array of dimension NCV.
  219. blas_int ldz = n;
  220. // seupd() will output directly into the eigval and eigvec objects.
  221. eigval.zeros(n_eigvals);
  222. eigvec.zeros(n, n_eigvals);
  223. arpack::seupd(&rvec, &howmny, select.memptr(), eigval.memptr(), eigvec.memptr(), &ldz, (eT*) NULL, &bmat, &n, which, &nev, &tol, resid.memptr(), &ncv, v.memptr(), &ldv, iparam.memptr(), ipntr.memptr(), workd.memptr(), workl.memptr(), &lworkl, &info);
  224. // Check for errors.
  225. if(info != 0)
  226. {
  227. arma_debug_warn("eigs_sym(): ARPACK error ", info, " in seupd()");
  228. return false;
  229. }
  230. return (info == 0);
  231. }
  232. #else
  233. {
  234. arma_ignore(eigval);
  235. arma_ignore(eigvec);
  236. arma_ignore(X);
  237. arma_ignore(n_eigvals);
  238. arma_ignore(form_str);
  239. arma_ignore(default_tol);
  240. return false;
  241. }
  242. #endif
  243. }
  244. //! immediate eigendecomposition of non-symmetric real sparse object
  245. template<typename T, typename T1>
  246. inline
  247. bool
  248. sp_auxlib::eigs_gen(Col< std::complex<T> >& eigval, Mat< std::complex<T> >& eigvec, const SpBase<T, T1>& X, const uword n_eigvals, const char* form_str, const T default_tol)
  249. {
  250. arma_extra_debug_sigprint();
  251. #if defined(ARMA_USE_NEWARP)
  252. {
  253. const unwrap_spmat<T1> U(X.get_ref());
  254. return sp_auxlib::eigs_gen_newarp(eigval, eigvec, U.M, n_eigvals, form_str, default_tol);
  255. }
  256. #elif defined(ARMA_USE_ARPACK)
  257. {
  258. const unwrap_spmat<T1> U(X.get_ref());
  259. return sp_auxlib::eigs_gen_arpack(eigval, eigvec, U.M, n_eigvals, form_str, default_tol);
  260. }
  261. #else
  262. {
  263. arma_ignore(eigval);
  264. arma_ignore(eigvec);
  265. arma_ignore(X);
  266. arma_ignore(n_eigvals);
  267. arma_ignore(form_str);
  268. arma_ignore(default_tol);
  269. arma_stop_logic_error("eigs_gen(): use of NEWARP or ARPACK must be enabled");
  270. return false;
  271. }
  272. #endif
  273. }
  274. template<typename T>
  275. inline
  276. bool
  277. sp_auxlib::eigs_gen_newarp(Col< std::complex<T> >& eigval, Mat< std::complex<T> >& eigvec, const SpMat<T>& X, const uword n_eigvals, const char* form_str, const T default_tol)
  278. {
  279. arma_extra_debug_sigprint();
  280. #if defined(ARMA_USE_NEWARP)
  281. {
  282. const form_type form_val = sp_auxlib::interpret_form_str(form_str);
  283. arma_debug_check( (form_val == form_none), "eigs_gen(): unknown form specified" );
  284. const newarp::SparseGenMatProd<T> op(X);
  285. arma_debug_check( (op.n_rows != op.n_cols), "eigs_gen(): given matrix must be square sized" );
  286. arma_debug_check( (n_eigvals + 1 >= op.n_rows), "eigs_gen(): n_eigvals + 1 must be less than the number of rows in the matrix" );
  287. // If the matrix is empty, the case is trivial.
  288. if( (op.n_cols == 0) || (n_eigvals == 0) ) // We already know n_cols == n_rows.
  289. {
  290. eigval.reset();
  291. eigvec.reset();
  292. return true;
  293. }
  294. uword n = op.n_rows;
  295. uword ncv = n_eigvals + 2 + 1;
  296. if(ncv < (2 * n_eigvals + 1)) { ncv = 2 * n_eigvals + 1; }
  297. if(ncv > n) { ncv = n; }
  298. T tol = (std::max)(default_tol, std::numeric_limits<T>::epsilon());
  299. // eigval.set_size(n_eigvals);
  300. // eigvec.set_size(n, n_eigvals);
  301. bool status = true;
  302. uword nconv = 0;
  303. try
  304. {
  305. if(form_val == form_lm)
  306. {
  307. newarp::GenEigsSolver< T, newarp::EigsSelect::LARGEST_MAGN, newarp::SparseGenMatProd<T> > eigs(op, n_eigvals, ncv);
  308. eigs.init();
  309. nconv = eigs.compute(1000, tol);
  310. eigval = eigs.eigenvalues();
  311. eigvec = eigs.eigenvectors();
  312. }
  313. else
  314. if(form_val == form_sm)
  315. {
  316. newarp::GenEigsSolver< T, newarp::EigsSelect::SMALLEST_MAGN, newarp::SparseGenMatProd<T> > eigs(op, n_eigvals, ncv);
  317. eigs.init();
  318. nconv = eigs.compute(1000, tol);
  319. eigval = eigs.eigenvalues();
  320. eigvec = eigs.eigenvectors();
  321. }
  322. else
  323. if(form_val == form_lr)
  324. {
  325. newarp::GenEigsSolver< T, newarp::EigsSelect::LARGEST_REAL, newarp::SparseGenMatProd<T> > eigs(op, n_eigvals, ncv);
  326. eigs.init();
  327. nconv = eigs.compute(1000, tol);
  328. eigval = eigs.eigenvalues();
  329. eigvec = eigs.eigenvectors();
  330. }
  331. else
  332. if(form_val == form_sr)
  333. {
  334. newarp::GenEigsSolver< T, newarp::EigsSelect::SMALLEST_REAL, newarp::SparseGenMatProd<T> > eigs(op, n_eigvals, ncv);
  335. eigs.init();
  336. nconv = eigs.compute(1000, tol);
  337. eigval = eigs.eigenvalues();
  338. eigvec = eigs.eigenvectors();
  339. }
  340. else
  341. if(form_val == form_li)
  342. {
  343. newarp::GenEigsSolver< T, newarp::EigsSelect::LARGEST_IMAG, newarp::SparseGenMatProd<T> > eigs(op, n_eigvals, ncv);
  344. eigs.init();
  345. nconv = eigs.compute(1000, tol);
  346. eigval = eigs.eigenvalues();
  347. eigvec = eigs.eigenvectors();
  348. }
  349. else
  350. if(form_val == form_si)
  351. {
  352. newarp::GenEigsSolver< T, newarp::EigsSelect::SMALLEST_IMAG, newarp::SparseGenMatProd<T> > eigs(op, n_eigvals, ncv);
  353. eigs.init();
  354. nconv = eigs.compute(1000, tol);
  355. eigval = eigs.eigenvalues();
  356. eigvec = eigs.eigenvectors();
  357. }
  358. }
  359. catch(const std::runtime_error&)
  360. {
  361. status = false;
  362. }
  363. if(status == true)
  364. {
  365. if(nconv == 0) { status = false; }
  366. }
  367. return status;
  368. }
  369. #else
  370. {
  371. arma_ignore(eigval);
  372. arma_ignore(eigvec);
  373. arma_ignore(X);
  374. arma_ignore(n_eigvals);
  375. arma_ignore(form_str);
  376. arma_ignore(default_tol);
  377. return false;
  378. }
  379. #endif
  380. }
  381. template<typename T>
  382. inline
  383. bool
  384. sp_auxlib::eigs_gen_arpack(Col< std::complex<T> >& eigval, Mat< std::complex<T> >& eigvec, const SpMat<T>& X, const uword n_eigvals, const char* form_str, const T default_tol)
  385. {
  386. arma_extra_debug_sigprint();
  387. #if defined(ARMA_USE_ARPACK)
  388. {
  389. const form_type form_val = sp_auxlib::interpret_form_str(form_str);
  390. arma_debug_check( (form_val == form_none), "eigs_gen(): unknown form specified" );
  391. char which_lm[3] = "LM";
  392. char which_sm[3] = "SM";
  393. char which_lr[3] = "LR";
  394. char which_sr[3] = "SR";
  395. char which_li[3] = "LI";
  396. char which_si[3] = "SI";
  397. char* which;
  398. switch(form_val)
  399. {
  400. case form_lm: which = which_lm; break;
  401. case form_sm: which = which_sm; break;
  402. case form_lr: which = which_lr; break;
  403. case form_sr: which = which_sr; break;
  404. case form_li: which = which_li; break;
  405. case form_si: which = which_si; break;
  406. default: which = which_lm;
  407. }
  408. // Make sure it's square.
  409. arma_debug_check( (X.n_rows != X.n_cols), "eigs_gen(): given matrix must be square sized" );
  410. // Make sure we aren't asking for every eigenvalue.
  411. arma_debug_check( (n_eigvals + 1 >= X.n_rows), "eigs_gen(): n_eigvals + 1 must be less than the number of rows in the matrix" );
  412. // If the matrix is empty, the case is trivial.
  413. if( (X.n_cols == 0) || (n_eigvals == 0) ) // We already know n_cols == n_rows.
  414. {
  415. eigval.reset();
  416. eigvec.reset();
  417. return true;
  418. }
  419. // Set up variables that get used for neupd().
  420. blas_int n, ncv, ldv, lworkl, info;
  421. T tol = default_tol;
  422. podarray<T> resid, v, workd, workl;
  423. podarray<blas_int> iparam, ipntr;
  424. podarray<T> rwork; // Not used in the real case.
  425. run_aupd(n_eigvals, which, X, false /* gen, not sym */, n, tol, resid, ncv, v, ldv, iparam, ipntr, workd, workl, lworkl, rwork, info);
  426. if(info != 0)
  427. {
  428. return false;
  429. }
  430. // The process has converged, and now we need to recover the actual eigenvectors using neupd().
  431. blas_int rvec = 1; // .TRUE
  432. blas_int nev = n_eigvals;
  433. char howmny = 'A';
  434. char bmat = 'I'; // We are considering the standard eigenvalue problem.
  435. podarray<blas_int> select(ncv); // Logical array of dimension NCV.
  436. podarray<T> dr(nev + 1); // Real array of dimension NEV + 1.
  437. podarray<T> di(nev + 1); // Real array of dimension NEV + 1.
  438. podarray<T> z(n * (nev + 1)); // Real N by NEV array if HOWMNY = 'A'.
  439. blas_int ldz = n;
  440. podarray<T> workev(3 * ncv);
  441. dr.zeros();
  442. di.zeros();
  443. z.zeros();
  444. arpack::neupd(&rvec, &howmny, select.memptr(), dr.memptr(), di.memptr(), z.memptr(), &ldz, (T*) NULL, (T*) NULL, workev.memptr(), &bmat, &n, which, &nev, &tol, resid.memptr(), &ncv, v.memptr(), &ldv, iparam.memptr(), ipntr.memptr(), workd.memptr(), workl.memptr(), &lworkl, rwork.memptr(), &info);
  445. // Check for errors.
  446. if(info != 0)
  447. {
  448. arma_debug_warn("eigs_gen(): ARPACK error ", info, " in neupd()");
  449. return false;
  450. }
  451. // Put it into the outputs.
  452. eigval.set_size(n_eigvals);
  453. eigvec.zeros(n, n_eigvals);
  454. for (uword i = 0; i < n_eigvals; ++i)
  455. {
  456. eigval[i] = std::complex<T>(dr[i], di[i]);
  457. }
  458. // Now recover the eigenvectors.
  459. for (uword i = 0; i < n_eigvals; ++i)
  460. {
  461. // ARPACK ?neupd lays things out kinda odd in memory;
  462. // so does LAPACK ?geev -- see auxlib::eig_gen()
  463. if((i < n_eigvals - 1) && (eigval[i] == std::conj(eigval[i + 1])))
  464. {
  465. for (uword j = 0; j < uword(n); ++j)
  466. {
  467. eigvec.at(j, i) = std::complex<T>(z[n * i + j], z[n * (i + 1) + j]);
  468. eigvec.at(j, i + 1) = std::complex<T>(z[n * i + j], -z[n * (i + 1) + j]);
  469. }
  470. ++i; // Skip the next one.
  471. }
  472. else
  473. if((i == n_eigvals - 1) && (std::complex<T>(eigval[i]).imag() != 0.0))
  474. {
  475. // We don't have the matched conjugate eigenvalue.
  476. for (uword j = 0; j < uword(n); ++j)
  477. {
  478. eigvec.at(j, i) = std::complex<T>(z[n * i + j], z[n * (i + 1) + j]);
  479. }
  480. }
  481. else
  482. {
  483. // The eigenvector is entirely real.
  484. for (uword j = 0; j < uword(n); ++j)
  485. {
  486. eigvec.at(j, i) = std::complex<T>(z[n * i + j], T(0));
  487. }
  488. }
  489. }
  490. return (info == 0);
  491. }
  492. #else
  493. {
  494. arma_ignore(eigval);
  495. arma_ignore(eigvec);
  496. arma_ignore(X);
  497. arma_ignore(n_eigvals);
  498. arma_ignore(form_str);
  499. arma_ignore(default_tol);
  500. return false;
  501. }
  502. #endif
  503. }
  504. //! immediate eigendecomposition of non-symmetric complex sparse object
  505. template<typename T, typename T1>
  506. inline
  507. bool
  508. sp_auxlib::eigs_gen(Col< std::complex<T> >& eigval, Mat< std::complex<T> >& eigvec, const SpBase< std::complex<T>, T1>& X_expr, const uword n_eigvals, const char* form_str, const T default_tol)
  509. {
  510. arma_extra_debug_sigprint();
  511. #if defined(ARMA_USE_ARPACK)
  512. {
  513. typedef typename std::complex<T> eT;
  514. const form_type form_val = sp_auxlib::interpret_form_str(form_str);
  515. arma_debug_check( (form_val == form_none), "eigs_gen(): unknown form specified" );
  516. char which_lm[3] = "LM";
  517. char which_sm[3] = "SM";
  518. char which_lr[3] = "LR";
  519. char which_sr[3] = "SR";
  520. char which_li[3] = "LI";
  521. char which_si[3] = "SI";
  522. char* which;
  523. switch(form_val)
  524. {
  525. case form_lm: which = which_lm; break;
  526. case form_sm: which = which_sm; break;
  527. case form_lr: which = which_lr; break;
  528. case form_sr: which = which_sr; break;
  529. case form_li: which = which_li; break;
  530. case form_si: which = which_si; break;
  531. default: which = which_lm;
  532. }
  533. const unwrap_spmat<T1> U(X_expr.get_ref());
  534. const SpMat<eT>& X = U.M;
  535. // Make sure it's square.
  536. arma_debug_check( (X.n_rows != X.n_cols), "eigs_gen(): given matrix must be square sized" );
  537. // Make sure we aren't asking for every eigenvalue.
  538. arma_debug_check( (n_eigvals + 1 >= X.n_rows), "eigs_gen(): n_eigvals + 1 must be less than the number of rows in the matrix" );
  539. // If the matrix is empty, the case is trivial.
  540. if( (X.n_cols == 0) || (n_eigvals == 0) ) // We already know n_cols == n_rows.
  541. {
  542. eigval.reset();
  543. eigvec.reset();
  544. return true;
  545. }
  546. // Set up variables that get used for neupd().
  547. blas_int n, ncv, ldv, lworkl, info;
  548. T tol = default_tol;
  549. podarray< std::complex<T> > resid, v, workd, workl;
  550. podarray<blas_int> iparam, ipntr;
  551. podarray<T> rwork;
  552. run_aupd(n_eigvals, which, X, false /* gen, not sym */, n, tol, resid, ncv, v, ldv, iparam, ipntr, workd, workl, lworkl, rwork, info);
  553. if(info != 0)
  554. {
  555. return false;
  556. }
  557. // The process has converged, and now we need to recover the actual eigenvectors using neupd().
  558. blas_int rvec = 1; // .TRUE
  559. blas_int nev = n_eigvals;
  560. char howmny = 'A';
  561. char bmat = 'I'; // We are considering the standard eigenvalue problem.
  562. podarray<blas_int> select(ncv); // Logical array of dimension NCV.
  563. podarray<std::complex<T> > d(nev + 1); // Real array of dimension NEV + 1.
  564. podarray<std::complex<T> > z(n * nev); // Real N by NEV array if HOWMNY = 'A'.
  565. blas_int ldz = n;
  566. podarray<std::complex<T> > workev(2 * ncv);
  567. // Prepare the outputs; neupd() will write directly to them.
  568. eigval.zeros(n_eigvals);
  569. eigvec.zeros(n, n_eigvals);
  570. std::complex<T> sigma;
  571. arpack::neupd(&rvec, &howmny, select.memptr(), eigval.memptr(),
  572. (std::complex<T>*) NULL, eigvec.memptr(), &ldz, (std::complex<T>*) &sigma, (std::complex<T>*) NULL, workev.memptr(), &bmat, &n, which, &nev, &tol, resid.memptr(), &ncv, v.memptr(), &ldv, iparam.memptr(), ipntr.memptr(), workd.memptr(), workl.memptr(), &lworkl, rwork.memptr(), &info);
  573. // Check for errors.
  574. if(info != 0)
  575. {
  576. arma_debug_warn("eigs_gen(): ARPACK error ", info, " in neupd()");
  577. return false;
  578. }
  579. return (info == 0);
  580. }
  581. #else
  582. {
  583. arma_ignore(eigval);
  584. arma_ignore(eigvec);
  585. arma_ignore(X_expr);
  586. arma_ignore(n_eigvals);
  587. arma_ignore(form_str);
  588. arma_ignore(default_tol);
  589. arma_stop_logic_error("eigs_gen(): use of ARPACK must be enabled for decomposition of complex matrices");
  590. return false;
  591. }
  592. #endif
  593. }
  594. template<typename T1, typename T2>
  595. inline
  596. bool
  597. sp_auxlib::spsolve_simple(Mat<typename T1::elem_type>& X, const SpBase<typename T1::elem_type, T1>& A_expr, const Base<typename T1::elem_type, T2>& B_expr, const superlu_opts& user_opts)
  598. {
  599. arma_extra_debug_sigprint();
  600. #if defined(ARMA_USE_SUPERLU)
  601. {
  602. typedef typename T1::elem_type eT;
  603. superlu::superlu_options_t options;
  604. sp_auxlib::set_superlu_opts(options, user_opts);
  605. const unwrap_spmat<T1> tmp1(A_expr.get_ref());
  606. const SpMat<eT>& A = tmp1.M;
  607. X = B_expr.get_ref(); // superlu::gssv() uses X as input (the B matrix) and as output (the solution)
  608. if(A.n_rows > A.n_cols)
  609. {
  610. arma_stop_logic_error("spsolve(): solving over-determined systems currently not supported");
  611. X.soft_reset();
  612. return false;
  613. }
  614. else if(A.n_rows < A.n_cols)
  615. {
  616. arma_stop_logic_error("spsolve(): solving under-determined systems currently not supported");
  617. X.soft_reset();
  618. return false;
  619. }
  620. arma_debug_check( (A.n_rows != X.n_rows), "spsolve(): number of rows in the given objects must be the same" );
  621. if(A.is_empty() || X.is_empty())
  622. {
  623. X.zeros(A.n_cols, X.n_cols);
  624. return true;
  625. }
  626. if(A.n_nonzero == uword(0)) { X.soft_reset(); return false; }
  627. if(arma_config::debug)
  628. {
  629. bool overflow;
  630. overflow = (A.n_nonzero > INT_MAX);
  631. overflow = (A.n_rows > INT_MAX) || overflow;
  632. overflow = (A.n_cols > INT_MAX) || overflow;
  633. overflow = (X.n_rows > INT_MAX) || overflow;
  634. overflow = (X.n_cols > INT_MAX) || overflow;
  635. if(overflow)
  636. {
  637. arma_stop_runtime_error("spsolve(): integer overflow: matrix dimensions are too large for integer type used by SuperLU");
  638. return false;
  639. }
  640. }
  641. superlu::SuperMatrix x; arrayops::inplace_set(reinterpret_cast<char*>(&x), char(0), sizeof(superlu::SuperMatrix));
  642. superlu::SuperMatrix a; arrayops::inplace_set(reinterpret_cast<char*>(&a), char(0), sizeof(superlu::SuperMatrix));
  643. const bool status_x = wrap_to_supermatrix(x, X);
  644. const bool status_a = copy_to_supermatrix(a, A);
  645. if( (status_x == false) || (status_a == false) )
  646. {
  647. destroy_supermatrix(a);
  648. destroy_supermatrix(x);
  649. X.soft_reset();
  650. return false;
  651. }
  652. superlu::SuperMatrix l; arrayops::inplace_set(reinterpret_cast<char*>(&l), char(0), sizeof(superlu::SuperMatrix));
  653. superlu::SuperMatrix u; arrayops::inplace_set(reinterpret_cast<char*>(&u), char(0), sizeof(superlu::SuperMatrix));
  654. // paranoia: use SuperLU's memory allocation, in case it reallocs
  655. int* perm_c = (int*) superlu::malloc( (A.n_cols+1) * sizeof(int)); // extra paranoia: increase array length by 1
  656. int* perm_r = (int*) superlu::malloc( (A.n_rows+1) * sizeof(int));
  657. arma_check_bad_alloc( (perm_c == 0), "spsolve(): out of memory" );
  658. arma_check_bad_alloc( (perm_r == 0), "spsolve(): out of memory" );
  659. arrayops::inplace_set(perm_c, 0, A.n_cols+1);
  660. arrayops::inplace_set(perm_r, 0, A.n_rows+1);
  661. superlu::SuperLUStat_t stat;
  662. superlu::init_stat(&stat);
  663. int info = 0; // Return code.
  664. arma_extra_debug_print("superlu::gssv()");
  665. superlu::gssv<eT>(&options, &a, perm_c, perm_r, &l, &u, &x, &stat, &info);
  666. // Process the return code.
  667. if( (info > 0) && (info <= int(A.n_cols)) )
  668. {
  669. // std::ostringstream tmp;
  670. // tmp << "spsolve(): could not solve system; LU factorisation completed, but detected zero in U(" << (info-1) << ',' << (info-1) << ')';
  671. // arma_debug_warn(tmp.str());
  672. }
  673. else
  674. if(info > int(A.n_cols))
  675. {
  676. arma_debug_warn("spsolve(): memory allocation failure: could not allocate ", (info - int(A.n_cols)), " bytes");
  677. }
  678. else
  679. if(info < 0)
  680. {
  681. arma_debug_warn("spsolve(): unknown SuperLU error code from gssv(): ", info);
  682. }
  683. superlu::free_stat(&stat);
  684. superlu::free(perm_c);
  685. superlu::free(perm_r);
  686. destroy_supermatrix(u);
  687. destroy_supermatrix(l);
  688. destroy_supermatrix(a);
  689. destroy_supermatrix(x); // No need to extract the data from x, since it's using the same memory as X
  690. return (info == 0);
  691. }
  692. #else
  693. {
  694. arma_ignore(X);
  695. arma_ignore(A_expr);
  696. arma_ignore(B_expr);
  697. arma_ignore(user_opts);
  698. arma_stop_logic_error("spsolve(): use of SuperLU must be enabled");
  699. return false;
  700. }
  701. #endif
  702. }
  703. template<typename T1, typename T2>
  704. inline
  705. bool
  706. sp_auxlib::spsolve_refine(Mat<typename T1::elem_type>& X, typename T1::pod_type& out_rcond, const SpBase<typename T1::elem_type, T1>& A_expr, const Base<typename T1::elem_type, T2>& B_expr, const superlu_opts& user_opts)
  707. {
  708. arma_extra_debug_sigprint();
  709. #if defined(ARMA_USE_SUPERLU)
  710. {
  711. typedef typename T1::pod_type T;
  712. typedef typename T1::elem_type eT;
  713. superlu::superlu_options_t options;
  714. sp_auxlib::set_superlu_opts(options, user_opts);
  715. const unwrap_spmat<T1> tmp1(A_expr.get_ref());
  716. const SpMat<eT>& A = tmp1.M;
  717. const unwrap<T2> tmp2(B_expr.get_ref());
  718. const Mat<eT>& B_unwrap = tmp2.M;
  719. const bool B_is_modified = ( (user_opts.equilibrate) || (&B_unwrap == &X) );
  720. Mat<eT> B_copy; if(B_is_modified) { B_copy = B_unwrap; }
  721. const Mat<eT>& B = (B_is_modified) ? B_copy : B_unwrap;
  722. if(A.n_rows > A.n_cols)
  723. {
  724. arma_stop_logic_error("spsolve(): solving over-determined systems currently not supported");
  725. X.soft_reset();
  726. return false;
  727. }
  728. else if(A.n_rows < A.n_cols)
  729. {
  730. arma_stop_logic_error("spsolve(): solving under-determined systems currently not supported");
  731. X.soft_reset();
  732. return false;
  733. }
  734. arma_debug_check( (A.n_rows != B.n_rows), "spsolve(): number of rows in the given objects must be the same" );
  735. X.zeros(A.n_cols, B.n_cols); // set the elements to zero, as we don't trust the SuperLU spaghetti code
  736. if(A.is_empty() || B.is_empty())
  737. {
  738. return true;
  739. }
  740. if(A.n_nonzero == uword(0)) { X.soft_reset(); return false; }
  741. if(arma_config::debug)
  742. {
  743. bool overflow;
  744. overflow = (A.n_nonzero > INT_MAX);
  745. overflow = (A.n_rows > INT_MAX) || overflow;
  746. overflow = (A.n_cols > INT_MAX) || overflow;
  747. overflow = (B.n_rows > INT_MAX) || overflow;
  748. overflow = (B.n_cols > INT_MAX) || overflow;
  749. overflow = (X.n_rows > INT_MAX) || overflow;
  750. overflow = (X.n_cols > INT_MAX) || overflow;
  751. if(overflow)
  752. {
  753. arma_stop_runtime_error("spsolve(): integer overflow: matrix dimensions are too large for integer type used by SuperLU");
  754. return false;
  755. }
  756. }
  757. superlu::SuperMatrix x; arrayops::inplace_set(reinterpret_cast<char*>(&x), char(0), sizeof(superlu::SuperMatrix));
  758. superlu::SuperMatrix a; arrayops::inplace_set(reinterpret_cast<char*>(&a), char(0), sizeof(superlu::SuperMatrix));
  759. superlu::SuperMatrix b; arrayops::inplace_set(reinterpret_cast<char*>(&b), char(0), sizeof(superlu::SuperMatrix));
  760. const bool status_x = wrap_to_supermatrix(x, X);
  761. const bool status_a = copy_to_supermatrix(a, A); // NOTE: superlu::gssvx() modifies 'a' if equilibration is enabled
  762. const bool status_b = wrap_to_supermatrix(b, B); // NOTE: superlu::gssvx() modifies 'b' if equilibration is enabled
  763. if( (status_x == false) || (status_a == false) || (status_b == false) )
  764. {
  765. destroy_supermatrix(x);
  766. destroy_supermatrix(a);
  767. destroy_supermatrix(b);
  768. X.soft_reset();
  769. return false;
  770. }
  771. superlu::SuperMatrix l; arrayops::inplace_set(reinterpret_cast<char*>(&l), char(0), sizeof(superlu::SuperMatrix));
  772. superlu::SuperMatrix u; arrayops::inplace_set(reinterpret_cast<char*>(&u), char(0), sizeof(superlu::SuperMatrix));
  773. // paranoia: use SuperLU's memory allocation, in case it reallocs
  774. int* perm_c = (int*) superlu::malloc( (A.n_cols+1) * sizeof(int) ); // extra paranoia: increase array length by 1
  775. int* perm_r = (int*) superlu::malloc( (A.n_rows+1) * sizeof(int) );
  776. int* etree = (int*) superlu::malloc( (A.n_cols+1) * sizeof(int) );
  777. T* R = (T*) superlu::malloc( (A.n_rows+1) * sizeof(T) );
  778. T* C = (T*) superlu::malloc( (A.n_cols+1) * sizeof(T) );
  779. T* ferr = (T*) superlu::malloc( (B.n_cols+1) * sizeof(T) );
  780. T* berr = (T*) superlu::malloc( (B.n_cols+1) * sizeof(T) );
  781. arma_check_bad_alloc( (perm_c == 0), "spsolve(): out of memory" );
  782. arma_check_bad_alloc( (perm_r == 0), "spsolve(): out of memory" );
  783. arma_check_bad_alloc( (etree == 0), "spsolve(): out of memory" );
  784. arma_check_bad_alloc( (R == 0), "spsolve(): out of memory" );
  785. arma_check_bad_alloc( (C == 0), "spsolve(): out of memory" );
  786. arma_check_bad_alloc( (ferr == 0), "spsolve(): out of memory" );
  787. arma_check_bad_alloc( (berr == 0), "spsolve(): out of memory" );
  788. arrayops::inplace_set(perm_c, int(0), A.n_cols+1);
  789. arrayops::inplace_set(perm_r, int(0), A.n_rows+1);
  790. arrayops::inplace_set(etree, int(0), A.n_cols+1);
  791. arrayops::inplace_set(R, T(0), A.n_rows+1);
  792. arrayops::inplace_set(C, T(0), A.n_cols+1);
  793. arrayops::inplace_set(ferr, T(0), B.n_cols+1);
  794. arrayops::inplace_set(berr, T(0), B.n_cols+1);
  795. superlu::GlobalLU_t glu;
  796. arrayops::inplace_set(reinterpret_cast<char*>(&glu), char(0), sizeof(superlu::GlobalLU_t));
  797. superlu::mem_usage_t mu;
  798. arrayops::inplace_set(reinterpret_cast<char*>(&mu), char(0), sizeof(superlu::mem_usage_t));
  799. superlu::SuperLUStat_t stat;
  800. superlu::init_stat(&stat);
  801. char equed[8]; // extra characters for paranoia
  802. T rpg = T(0);
  803. T rcond = T(0);
  804. int info = int(0); // Return code.
  805. char work[8];
  806. int lwork = int(0); // 0 means superlu will allocate memory
  807. arma_extra_debug_print("superlu::gssvx()");
  808. superlu::gssvx<eT>(&options, &a, perm_c, perm_r, etree, equed, R, C, &l, &u, &work[0], lwork, &b, &x, &rpg, &rcond, ferr, berr, &glu, &mu, &stat, &info);
  809. bool status = false;
  810. // Process the return code.
  811. if(info == 0)
  812. {
  813. status = true;
  814. }
  815. if( (info > 0) && (info <= int(A.n_cols)) )
  816. {
  817. // std::ostringstream tmp;
  818. // tmp << "spsolve(): could not solve system; LU factorisation completed, but detected zero in U(" << (info-1) << ',' << (info-1) << ')';
  819. // arma_debug_warn(tmp.str());
  820. }
  821. else
  822. if( (info == int(A.n_cols+1)) && (user_opts.allow_ugly) )
  823. {
  824. arma_debug_warn("spsolve(): system is singular to working precision (rcond: ", rcond, ")");
  825. status = true;
  826. }
  827. else
  828. if(info > int(A.n_cols+1))
  829. {
  830. arma_debug_warn("spsolve(): memory allocation failure: could not allocate ", (info - int(A.n_cols)), " bytes");
  831. }
  832. else
  833. if(info < 0)
  834. {
  835. arma_debug_warn("spsolve(): unknown SuperLU error code from gssvx(): ", info);
  836. }
  837. superlu::free_stat(&stat);
  838. superlu::free(berr);
  839. superlu::free(ferr);
  840. superlu::free(C);
  841. superlu::free(R);
  842. superlu::free(etree);
  843. superlu::free(perm_r);
  844. superlu::free(perm_c);
  845. destroy_supermatrix(u);
  846. destroy_supermatrix(l);
  847. destroy_supermatrix(b);
  848. destroy_supermatrix(a);
  849. destroy_supermatrix(x); // No need to extract the data from x, since it's using the same memory as X
  850. out_rcond = rcond;
  851. return status;
  852. }
  853. #else
  854. {
  855. arma_ignore(X);
  856. arma_ignore(out_rcond);
  857. arma_ignore(A_expr);
  858. arma_ignore(B_expr);
  859. arma_ignore(user_opts);
  860. arma_stop_logic_error("spsolve(): use of SuperLU must be enabled");
  861. return false;
  862. }
  863. #endif
  864. }
  865. #if defined(ARMA_USE_SUPERLU)
  866. inline
  867. void
  868. sp_auxlib::set_superlu_opts(superlu::superlu_options_t& options, const superlu_opts& user_opts)
  869. {
  870. arma_extra_debug_sigprint();
  871. // default options as the starting point
  872. superlu::set_default_opts(&options);
  873. // our settings
  874. options.Trans = superlu::NOTRANS;
  875. options.ConditionNumber = superlu::YES;
  876. // process user_opts
  877. if(user_opts.equilibrate == true) { options.Equil = superlu::YES; }
  878. if(user_opts.equilibrate == false) { options.Equil = superlu::NO; }
  879. if(user_opts.symmetric == true) { options.SymmetricMode = superlu::YES; }
  880. if(user_opts.symmetric == false) { options.SymmetricMode = superlu::NO; }
  881. options.DiagPivotThresh = user_opts.pivot_thresh;
  882. if(user_opts.permutation == superlu_opts::NATURAL) { options.ColPerm = superlu::NATURAL; }
  883. if(user_opts.permutation == superlu_opts::MMD_ATA) { options.ColPerm = superlu::MMD_ATA; }
  884. if(user_opts.permutation == superlu_opts::MMD_AT_PLUS_A) { options.ColPerm = superlu::MMD_AT_PLUS_A; }
  885. if(user_opts.permutation == superlu_opts::COLAMD) { options.ColPerm = superlu::COLAMD; }
  886. if(user_opts.refine == superlu_opts::REF_NONE) { options.IterRefine = superlu::NOREFINE; }
  887. if(user_opts.refine == superlu_opts::REF_SINGLE) { options.IterRefine = superlu::SLU_SINGLE; }
  888. if(user_opts.refine == superlu_opts::REF_DOUBLE) { options.IterRefine = superlu::SLU_DOUBLE; }
  889. if(user_opts.refine == superlu_opts::REF_EXTRA) { options.IterRefine = superlu::SLU_EXTRA; }
  890. }
  891. template<typename eT>
  892. inline
  893. bool
  894. sp_auxlib::copy_to_supermatrix(superlu::SuperMatrix& out, const SpMat<eT>& A)
  895. {
  896. arma_extra_debug_sigprint();
  897. // We store in column-major CSC.
  898. out.Stype = superlu::SLU_NC;
  899. if(is_float<eT>::value)
  900. {
  901. out.Dtype = superlu::SLU_S;
  902. }
  903. else
  904. if(is_double<eT>::value)
  905. {
  906. out.Dtype = superlu::SLU_D;
  907. }
  908. else
  909. if(is_cx_float<eT>::value)
  910. {
  911. out.Dtype = superlu::SLU_C;
  912. }
  913. else
  914. if(is_cx_double<eT>::value)
  915. {
  916. out.Dtype = superlu::SLU_Z;
  917. }
  918. out.Mtype = superlu::SLU_GE; // Just a general matrix. We don't know more now.
  919. // We have to actually create the object which stores the data.
  920. // This gets cleaned by destroy_supermatrix().
  921. // We have to use SuperLU's problematic memory allocation routines since they are
  922. // not guaranteed to be new and delete. See the comments in def_superlu.hpp
  923. superlu::NCformat* nc = (superlu::NCformat*)superlu::malloc(sizeof(superlu::NCformat));
  924. if(nc == NULL) { return false; }
  925. A.sync();
  926. nc->nnz = A.n_nonzero;
  927. nc->nzval = (void*) superlu::malloc(sizeof(eT) * A.n_nonzero );
  928. nc->colptr = (superlu::int_t*)superlu::malloc(sizeof(superlu::int_t) * (A.n_cols + 1));
  929. nc->rowind = (superlu::int_t*)superlu::malloc(sizeof(superlu::int_t) * A.n_nonzero );
  930. if( (nc->nzval == NULL) || (nc->colptr == NULL) || (nc->rowind == NULL) ) { return false; }
  931. // Fill the matrix.
  932. arrayops::copy((eT*) nc->nzval, A.values, A.n_nonzero);
  933. // // These have to be copied by hand, because the types may differ.
  934. // for (uword i = 0; i <= A.n_cols; ++i) { nc->colptr[i] = (int_t) A.col_ptrs[i]; }
  935. // for (uword i = 0; i < A.n_nonzero; ++i) { nc->rowind[i] = (int_t) A.row_indices[i]; }
  936. arrayops::convert(nc->colptr, A.col_ptrs, A.n_cols+1 );
  937. arrayops::convert(nc->rowind, A.row_indices, A.n_nonzero);
  938. out.nrow = A.n_rows;
  939. out.ncol = A.n_cols;
  940. out.Store = (void*) nc;
  941. return true;
  942. }
  943. template<typename eT>
  944. inline
  945. bool
  946. sp_auxlib::wrap_to_supermatrix(superlu::SuperMatrix& out, const Mat<eT>& A)
  947. {
  948. arma_extra_debug_sigprint();
  949. // NOTE: this function re-uses memory from matrix A
  950. // This is being stored as a dense matrix.
  951. out.Stype = superlu::SLU_DN;
  952. if(is_float<eT>::value)
  953. {
  954. out.Dtype = superlu::SLU_S;
  955. }
  956. else
  957. if(is_double<eT>::value)
  958. {
  959. out.Dtype = superlu::SLU_D;
  960. }
  961. else
  962. if(is_cx_float<eT>::value)
  963. {
  964. out.Dtype = superlu::SLU_C;
  965. }
  966. else
  967. if(is_cx_double<eT>::value)
  968. {
  969. out.Dtype = superlu::SLU_Z;
  970. }
  971. out.Mtype = superlu::SLU_GE;
  972. // We have to create the object that stores the data.
  973. superlu::DNformat* dn = (superlu::DNformat*)superlu::malloc(sizeof(superlu::DNformat));
  974. if(dn == NULL) { return false; }
  975. dn->lda = A.n_rows;
  976. dn->nzval = (void*) A.memptr(); // re-use memory instead of copying
  977. out.nrow = A.n_rows;
  978. out.ncol = A.n_cols;
  979. out.Store = (void*) dn;
  980. return true;
  981. }
  982. inline
  983. void
  984. sp_auxlib::destroy_supermatrix(superlu::SuperMatrix& out)
  985. {
  986. arma_extra_debug_sigprint();
  987. // Clean up.
  988. if(out.Stype == superlu::SLU_NC)
  989. {
  990. superlu::destroy_compcol_mat(&out);
  991. }
  992. else
  993. if(out.Stype == superlu::SLU_DN)
  994. {
  995. // superlu::destroy_dense_mat(&out);
  996. // since dn->nzval is set to re-use memory from a Mat object (which manages its own memory),
  997. // we cannot simply call superlu::destroy_dense_mat().
  998. // Only the out.Store structure can be freed.
  999. superlu::DNformat* dn = (superlu::DNformat*) out.Store;
  1000. if(dn != NULL) { superlu::free(dn); }
  1001. }
  1002. else
  1003. if(out.Stype == superlu::SLU_SC)
  1004. {
  1005. superlu::destroy_supernode_mat(&out);
  1006. }
  1007. else
  1008. {
  1009. // Uh, crap.
  1010. std::ostringstream tmp;
  1011. tmp << "sp_auxlib::destroy_supermatrix(): unhandled Stype" << std::endl;
  1012. tmp << "Stype val: " << out.Stype << std::endl;
  1013. tmp << "Stype name: ";
  1014. if(out.Stype == superlu::SLU_NC) { tmp << "SLU_NC"; }
  1015. if(out.Stype == superlu::SLU_NCP) { tmp << "SLU_NCP"; }
  1016. if(out.Stype == superlu::SLU_NR) { tmp << "SLU_NR"; }
  1017. if(out.Stype == superlu::SLU_SC) { tmp << "SLU_SC"; }
  1018. if(out.Stype == superlu::SLU_SCP) { tmp << "SLU_SCP"; }
  1019. if(out.Stype == superlu::SLU_SR) { tmp << "SLU_SR"; }
  1020. if(out.Stype == superlu::SLU_DN) { tmp << "SLU_DN"; }
  1021. if(out.Stype == superlu::SLU_NR_loc) { tmp << "SLU_NR_loc"; }
  1022. arma_debug_warn(tmp.str());
  1023. arma_stop_runtime_error("internal error: sp_auxlib::destroy_supermatrix()");
  1024. }
  1025. }
  1026. #endif
  1027. template<typename eT, typename T>
  1028. inline
  1029. void
  1030. sp_auxlib::run_aupd
  1031. (
  1032. const uword n_eigvals, char* which, const SpMat<T>& X, const bool sym,
  1033. blas_int& n, eT& tol,
  1034. podarray<T>& resid, blas_int& ncv, podarray<T>& v, blas_int& ldv,
  1035. podarray<blas_int>& iparam, podarray<blas_int>& ipntr,
  1036. podarray<T>& workd, podarray<T>& workl, blas_int& lworkl, podarray<eT>& rwork,
  1037. blas_int& info
  1038. )
  1039. {
  1040. #if defined(ARMA_USE_ARPACK)
  1041. {
  1042. // ARPACK provides a "reverse communication interface" which is an
  1043. // entertainingly archaic FORTRAN software engineering technique that
  1044. // basically means that we call saupd()/naupd() and it tells us with some
  1045. // return code what we need to do next (usually a matrix-vector product) and
  1046. // then call it again. So this results in some type of iterative process
  1047. // where we call saupd()/naupd() many times.
  1048. blas_int ido = 0; // This must be 0 for the first call.
  1049. char bmat = 'I'; // We are considering the standard eigenvalue problem.
  1050. n = X.n_rows; // The size of the matrix.
  1051. blas_int nev = n_eigvals;
  1052. resid.set_size(n);
  1053. // Two contraints on NCV: (NCV > NEV + 2) and (NCV <= N)
  1054. //
  1055. // We're calling either arpack::saupd() or arpack::naupd(),
  1056. // which have slighly different minimum constraint and recommended value for NCV:
  1057. // http://www.caam.rice.edu/software/ARPACK/UG/node136.html
  1058. // http://www.caam.rice.edu/software/ARPACK/UG/node138.html
  1059. ncv = nev + 2 + 1;
  1060. if (ncv < (2 * nev + 1)) { ncv = 2 * nev + 1; }
  1061. if (ncv > n ) { ncv = n; }
  1062. v.set_size(n * ncv); // Array N by NCV (output).
  1063. rwork.set_size(ncv); // Work array of size NCV for complex calls.
  1064. ldv = n; // "Leading dimension of V exactly as declared in the calling program."
  1065. // IPARAM: integer array of length 11.
  1066. iparam.zeros(11);
  1067. iparam(0) = 1; // Exact shifts (not provided by us).
  1068. iparam(2) = 1000; // Maximum iterations; all the examples use 300, but they were written in the ancient times.
  1069. iparam(6) = 1; // Mode 1: A * x = lambda * x.
  1070. // IPNTR: integer array of length 14 (output).
  1071. ipntr.set_size(14);
  1072. // Real work array used in the basic Arnoldi iteration for reverse communication.
  1073. workd.set_size(3 * n);
  1074. // lworkl must be at least 3 * NCV^2 + 6 * NCV.
  1075. lworkl = 3 * (ncv * ncv) + 6 * ncv;
  1076. // Real work array of length lworkl.
  1077. workl.set_size(lworkl);
  1078. info = 0; // Set to 0 initially to use random initial vector.
  1079. // All the parameters have been set or created. Time to loop a lot.
  1080. while (ido != 99)
  1081. {
  1082. // Call saupd() or naupd() with the current parameters.
  1083. if(sym)
  1084. arpack::saupd(&ido, &bmat, &n, which, &nev, &tol, resid.memptr(), &ncv, v.memptr(), &ldv, iparam.memptr(), ipntr.memptr(), workd.memptr(), workl.memptr(), &lworkl, &info);
  1085. else
  1086. arpack::naupd(&ido, &bmat, &n, which, &nev, &tol, resid.memptr(), &ncv, v.memptr(), &ldv, iparam.memptr(), ipntr.memptr(), workd.memptr(), workl.memptr(), &lworkl, rwork.memptr(), &info);
  1087. // What do we do now?
  1088. switch (ido)
  1089. {
  1090. case -1:
  1091. // fallthrough
  1092. case 1:
  1093. {
  1094. // We need to calculate the matrix-vector multiplication y = OP * x
  1095. // where x is of length n and starts at workd(ipntr(0)), and y is of
  1096. // length n and starts at workd(ipntr(1)).
  1097. // operator*(sp_mat, vec) doesn't properly put the result into the
  1098. // right place so we'll just reimplement it here for now...
  1099. // Set the output to point at the right memory. We have to subtract
  1100. // one from FORTRAN pointers...
  1101. Col<T> out(workd.memptr() + ipntr(1) - 1, n, false /* don't copy */);
  1102. // Set the input to point at the right memory.
  1103. Col<T> in(workd.memptr() + ipntr(0) - 1, n, false /* don't copy */);
  1104. out.zeros();
  1105. typename SpMat<T>::const_iterator x_it = X.begin();
  1106. typename SpMat<T>::const_iterator x_it_end = X.end();
  1107. while(x_it != x_it_end)
  1108. {
  1109. out[x_it.row()] += (*x_it) * in[x_it.col()];
  1110. ++x_it;
  1111. }
  1112. // No need to modify memory further since it was all done in-place.
  1113. break;
  1114. }
  1115. case 99:
  1116. // Nothing to do here, things have converged.
  1117. break;
  1118. default:
  1119. {
  1120. return; // Parent frame can look at the value of info.
  1121. }
  1122. }
  1123. }
  1124. // The process has ended; check the return code.
  1125. if( (info != 0) && (info != 1) )
  1126. {
  1127. // Print warnings if there was a failure.
  1128. if(sym)
  1129. {
  1130. arma_debug_warn("eigs_sym(): ARPACK error ", info, " in saupd()");
  1131. }
  1132. else
  1133. {
  1134. arma_debug_warn("eigs_gen(): ARPACK error ", info, " in naupd()");
  1135. }
  1136. return; // Parent frame can look at the value of info.
  1137. }
  1138. }
  1139. #else
  1140. arma_ignore(n_eigvals);
  1141. arma_ignore(which);
  1142. arma_ignore(X);
  1143. arma_ignore(sym);
  1144. arma_ignore(n);
  1145. arma_ignore(tol);
  1146. arma_ignore(resid);
  1147. arma_ignore(ncv);
  1148. arma_ignore(v);
  1149. arma_ignore(ldv);
  1150. arma_ignore(iparam);
  1151. arma_ignore(ipntr);
  1152. arma_ignore(workd);
  1153. arma_ignore(workl);
  1154. arma_ignore(lworkl);
  1155. arma_ignore(rwork);
  1156. arma_ignore(info);
  1157. #endif
  1158. }
  1159. template<typename eT>
  1160. inline
  1161. bool
  1162. sp_auxlib::rudimentary_sym_check(const SpMat<eT>& X)
  1163. {
  1164. arma_extra_debug_sigprint();
  1165. if(X.n_rows != X.n_cols) { return false; }
  1166. const eT tol = eT(10000) * std::numeric_limits<eT>::epsilon(); // allow some leeway
  1167. typename SpMat<eT>::const_iterator it = X.begin();
  1168. typename SpMat<eT>::const_iterator it_end = X.end();
  1169. const uword n_check_limit = (std::max)( uword(2), uword(X.n_nonzero/100) );
  1170. uword n_check = 1;
  1171. while( (it != it_end) && (n_check <= n_check_limit) )
  1172. {
  1173. const uword it_row = it.row();
  1174. const uword it_col = it.col();
  1175. if(it_row != it_col)
  1176. {
  1177. const eT A = (*it);
  1178. const eT B = X.at( it_col, it_row ); // deliberately swapped
  1179. const eT C = (std::max)(std::abs(A), std::abs(B));
  1180. const eT delta = std::abs(A - B);
  1181. if(( (delta <= tol) || (delta <= (C * tol)) ) == false) { return false; }
  1182. ++n_check;
  1183. }
  1184. ++it;
  1185. }
  1186. return true;
  1187. }
  1188. template<typename T>
  1189. inline
  1190. bool
  1191. sp_auxlib::rudimentary_sym_check(const SpMat< std::complex<T> >& X)
  1192. {
  1193. arma_extra_debug_sigprint();
  1194. // NOTE: the function name is a misnomer, as it checks for hermitian complex matrices;
  1195. // NOTE: for simplicity of use, the function name is the same as for real matrices
  1196. typedef typename std::complex<T> eT;
  1197. if(X.n_rows != X.n_cols) { return false; }
  1198. const T tol = T(10000) * std::numeric_limits<T>::epsilon(); // allow some leeway
  1199. typename SpMat<eT>::const_iterator it = X.begin();
  1200. typename SpMat<eT>::const_iterator it_end = X.end();
  1201. const uword n_check_limit = (std::max)( uword(2), uword(X.n_nonzero/100) );
  1202. uword n_check = 1;
  1203. while( (it != it_end) && (n_check <= n_check_limit) )
  1204. {
  1205. const uword it_row = it.row();
  1206. const uword it_col = it.col();
  1207. if(it_row != it_col)
  1208. {
  1209. const eT A = (*it);
  1210. const eT B = X.at( it_col, it_row ); // deliberately swapped
  1211. const T C_real = (std::max)(std::abs(A.real()), std::abs(B.real()));
  1212. const T C_imag = (std::max)(std::abs(A.imag()), std::abs(B.imag()));
  1213. const T delta_real = std::abs(A.real() - B.real());
  1214. const T delta_imag = std::abs(A.imag() + B.imag()); // take into account the conjugate
  1215. const bool okay_real = ( (delta_real <= tol) || (delta_real <= (C_real * tol)) );
  1216. const bool okay_imag = ( (delta_imag <= tol) || (delta_imag <= (C_imag * tol)) );
  1217. if( (okay_real == false) || (okay_imag == false) ) { return false; }
  1218. ++n_check;
  1219. }
  1220. else
  1221. {
  1222. const eT A = (*it);
  1223. if(std::abs(A.imag()) > tol) { return false; }
  1224. }
  1225. ++it;
  1226. }
  1227. return true;
  1228. }
  1229. //! @}