123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554 |
- // Copyright 2008-2016 Conrad Sanderson (http://conradsanderson.id.au)
- // Copyright 2008-2016 National ICT Australia (NICTA)
- //
- // Licensed under the Apache License, Version 2.0 (the "License");
- // you may not use this file except in compliance with the License.
- // You may obtain a copy of the License at
- // http://www.apache.org/licenses/LICENSE-2.0
- //
- // Unless required by applicable law or agreed to in writing, software
- // distributed under the License is distributed on an "AS IS" BASIS,
- // WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
- // See the License for the specific language governing permissions and
- // limitations under the License.
- // ------------------------------------------------------------------------
- //! \addtogroup sp_auxlib
- //! @{
- inline
- sp_auxlib::form_type
- sp_auxlib::interpret_form_str(const char* form_str)
- {
- arma_extra_debug_sigprint();
-
- // the order of the 3 if statements below is important
- if( form_str == NULL ) { return form_none; }
- if( form_str[0] == char(0) ) { return form_none; }
- if( form_str[1] == char(0) ) { return form_none; }
-
- const char c1 = form_str[0];
- const char c2 = form_str[1];
-
- if(c1 == 'l')
- {
- if(c2 == 'm') { return form_lm; }
- if(c2 == 'r') { return form_lr; }
- if(c2 == 'i') { return form_li; }
- if(c2 == 'a') { return form_la; }
- }
- else
- if(c1 == 's')
- {
- if(c2 == 'm') { return form_sm; }
- if(c2 == 'r') { return form_sr; }
- if(c2 == 'i') { return form_si; }
- if(c2 == 'a') { return form_sa; }
- }
-
- return form_none;
- }
- //! immediate eigendecomposition of symmetric real sparse object
- template<typename eT, typename T1>
- inline
- bool
- 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)
- {
- arma_extra_debug_sigprint();
-
- const unwrap_spmat<T1> U(X.get_ref());
-
- if((arma_config::debug) && (sp_auxlib::rudimentary_sym_check(U.M) == false))
- {
- if(is_cx<eT>::no ) { arma_debug_warn("eigs_sym(): given matrix is not symmetric"); }
- if(is_cx<eT>::yes) { arma_debug_warn("eigs_sym(): given matrix is not hermitian"); }
- }
-
- #if defined(ARMA_USE_NEWARP)
- {
- return sp_auxlib::eigs_sym_newarp(eigval, eigvec, U.M, n_eigvals, form_str, default_tol);
- }
- #elif defined(ARMA_USE_ARPACK)
- {
- return sp_auxlib::eigs_sym_arpack(eigval, eigvec, U.M, n_eigvals, form_str, default_tol);
- }
- #else
- {
- arma_ignore(eigval);
- arma_ignore(eigvec);
- arma_ignore(n_eigvals);
- arma_ignore(form_str);
- arma_ignore(default_tol);
-
- arma_stop_logic_error("eigs_sym(): use of NEWARP or ARPACK must be enabled");
- return false;
- }
- #endif
- }
- template<typename eT>
- inline
- bool
- 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)
- {
- arma_extra_debug_sigprint();
-
- #if defined(ARMA_USE_NEWARP)
- {
- const form_type form_val = sp_auxlib::interpret_form_str(form_str);
-
- arma_debug_check( (form_val != form_lm) && (form_val != form_sm) && (form_val != form_la) && (form_val != form_sa), "eigs_sym(): unknown form specified" );
-
- const newarp::SparseGenMatProd<eT> op(X);
-
- arma_debug_check( (op.n_rows != op.n_cols), "eigs_sym(): given matrix must be square sized" );
-
- arma_debug_check( (n_eigvals >= op.n_rows), "eigs_sym(): n_eigvals must be less than the number of rows in the matrix" );
-
- // If the matrix is empty, the case is trivial.
- if( (op.n_cols == 0) || (n_eigvals == 0) ) // We already know n_cols == n_rows.
- {
- eigval.reset();
- eigvec.reset();
- return true;
- }
-
- uword n = op.n_rows;
- uword ncv = n_eigvals + 2 + 1;
-
- if(ncv < (2 * n_eigvals + 1)) { ncv = 2 * n_eigvals + 1; }
- if(ncv > n) { ncv = n; }
-
- eT tol = (std::max)(default_tol, std::numeric_limits<eT>::epsilon());
-
- // eigval.set_size(n_eigvals);
- // eigvec.set_size(n, n_eigvals);
-
- bool status = true;
-
- uword nconv = 0;
-
- try
- {
- if(form_val == form_lm)
- {
- newarp::SymEigsSolver< eT, newarp::EigsSelect::LARGEST_MAGN, newarp::SparseGenMatProd<eT> > eigs(op, n_eigvals, ncv);
- eigs.init();
- nconv = eigs.compute(1000, tol);
- eigval = eigs.eigenvalues();
- eigvec = eigs.eigenvectors();
- }
- else
- if(form_val == form_sm)
- {
- newarp::SymEigsSolver< eT, newarp::EigsSelect::SMALLEST_MAGN, newarp::SparseGenMatProd<eT> > eigs(op, n_eigvals, ncv);
- eigs.init();
- nconv = eigs.compute(1000, tol);
- eigval = eigs.eigenvalues();
- eigvec = eigs.eigenvectors();
- }
- else
- if(form_val == form_la)
- {
- newarp::SymEigsSolver< eT, newarp::EigsSelect::LARGEST_ALGE, newarp::SparseGenMatProd<eT> > eigs(op, n_eigvals, ncv);
- eigs.init();
- nconv = eigs.compute(1000, tol);
- eigval = eigs.eigenvalues();
- eigvec = eigs.eigenvectors();
- }
- else
- if(form_val == form_sa)
- {
- newarp::SymEigsSolver< eT, newarp::EigsSelect::SMALLEST_ALGE, newarp::SparseGenMatProd<eT> > eigs(op, n_eigvals, ncv);
- eigs.init();
- nconv = eigs.compute(1000, tol);
- eigval = eigs.eigenvalues();
- eigvec = eigs.eigenvectors();
- }
- }
- catch(const std::runtime_error&)
- {
- status = false;
- }
-
- if(status == true)
- {
- if(nconv == 0) { status = false; }
- }
-
- return status;
- }
- #else
- {
- arma_ignore(eigval);
- arma_ignore(eigvec);
- arma_ignore(X);
- arma_ignore(n_eigvals);
- arma_ignore(form_str);
- arma_ignore(default_tol);
-
- return false;
- }
- #endif
- }
- template<typename eT>
- inline
- bool
- 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)
- {
- arma_extra_debug_sigprint();
-
- #if defined(ARMA_USE_ARPACK)
- {
- const form_type form_val = sp_auxlib::interpret_form_str(form_str);
-
- arma_debug_check( (form_val != form_lm) && (form_val != form_sm) && (form_val != form_la) && (form_val != form_sa), "eigs_sym(): unknown form specified" );
-
- char which_sm[3] = "SM";
- char which_lm[3] = "LM";
- char which_sa[3] = "SA";
- char which_la[3] = "LA";
- char* which;
- switch (form_val)
- {
- case form_sm: which = which_sm; break;
- case form_lm: which = which_lm; break;
- case form_sa: which = which_sa; break;
- case form_la: which = which_la; break;
- default: which = which_lm; break;
- }
-
- // Make sure it's square.
- arma_debug_check( (X.n_rows != X.n_cols), "eigs_sym(): given matrix must be square sized" );
-
- // Make sure we aren't asking for every eigenvalue.
- // The _saupd() functions allow asking for one more eigenvalue than the _naupd() functions.
- arma_debug_check( (n_eigvals >= X.n_rows), "eigs_sym(): n_eigvals must be less than the number of rows in the matrix" );
-
- // If the matrix is empty, the case is trivial.
- if( (X.n_cols == 0) || (n_eigvals == 0) ) // We already know n_cols == n_rows.
- {
- eigval.reset();
- eigvec.reset();
- return true;
- }
-
- // Set up variables that get used for neupd().
- blas_int n, ncv, ldv, lworkl, info;
- eT tol = default_tol;
- podarray<eT> resid, v, workd, workl;
- podarray<blas_int> iparam, ipntr;
- podarray<eT> rwork; // Not used in this case.
-
- run_aupd(n_eigvals, which, X, true /* sym, not gen */, n, tol, resid, ncv, v, ldv, iparam, ipntr, workd, workl, lworkl, rwork, info);
-
- if(info != 0)
- {
- return false;
- }
-
- // The process has converged, and now we need to recover the actual eigenvectors using seupd()
- blas_int rvec = 1; // .TRUE
- blas_int nev = n_eigvals;
-
- char howmny = 'A';
- char bmat = 'I'; // We are considering the standard eigenvalue problem.
-
- podarray<blas_int> select(ncv); // Logical array of dimension NCV.
- blas_int ldz = n;
-
- // seupd() will output directly into the eigval and eigvec objects.
- eigval.zeros(n_eigvals);
- eigvec.zeros(n, n_eigvals);
-
- 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);
-
- // Check for errors.
- if(info != 0)
- {
- arma_debug_warn("eigs_sym(): ARPACK error ", info, " in seupd()");
- return false;
- }
-
- return (info == 0);
- }
- #else
- {
- arma_ignore(eigval);
- arma_ignore(eigvec);
- arma_ignore(X);
- arma_ignore(n_eigvals);
- arma_ignore(form_str);
- arma_ignore(default_tol);
-
- return false;
- }
- #endif
- }
- //! immediate eigendecomposition of non-symmetric real sparse object
- template<typename T, typename T1>
- inline
- bool
- 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)
- {
- arma_extra_debug_sigprint();
-
- #if defined(ARMA_USE_NEWARP)
- {
- const unwrap_spmat<T1> U(X.get_ref());
-
- return sp_auxlib::eigs_gen_newarp(eigval, eigvec, U.M, n_eigvals, form_str, default_tol);
- }
- #elif defined(ARMA_USE_ARPACK)
- {
- const unwrap_spmat<T1> U(X.get_ref());
-
- return sp_auxlib::eigs_gen_arpack(eigval, eigvec, U.M, n_eigvals, form_str, default_tol);
- }
- #else
- {
- arma_ignore(eigval);
- arma_ignore(eigvec);
- arma_ignore(X);
- arma_ignore(n_eigvals);
- arma_ignore(form_str);
- arma_ignore(default_tol);
-
- arma_stop_logic_error("eigs_gen(): use of NEWARP or ARPACK must be enabled");
- return false;
- }
- #endif
- }
- template<typename T>
- inline
- bool
- 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)
- {
- arma_extra_debug_sigprint();
-
- #if defined(ARMA_USE_NEWARP)
- {
- const form_type form_val = sp_auxlib::interpret_form_str(form_str);
-
- arma_debug_check( (form_val == form_none), "eigs_gen(): unknown form specified" );
-
- const newarp::SparseGenMatProd<T> op(X);
-
- arma_debug_check( (op.n_rows != op.n_cols), "eigs_gen(): given matrix must be square sized" );
-
- 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" );
-
- // If the matrix is empty, the case is trivial.
- if( (op.n_cols == 0) || (n_eigvals == 0) ) // We already know n_cols == n_rows.
- {
- eigval.reset();
- eigvec.reset();
- return true;
- }
-
- uword n = op.n_rows;
- uword ncv = n_eigvals + 2 + 1;
-
- if(ncv < (2 * n_eigvals + 1)) { ncv = 2 * n_eigvals + 1; }
- if(ncv > n) { ncv = n; }
-
- T tol = (std::max)(default_tol, std::numeric_limits<T>::epsilon());
-
- // eigval.set_size(n_eigvals);
- // eigvec.set_size(n, n_eigvals);
-
- bool status = true;
-
- uword nconv = 0;
-
- try
- {
- if(form_val == form_lm)
- {
- newarp::GenEigsSolver< T, newarp::EigsSelect::LARGEST_MAGN, newarp::SparseGenMatProd<T> > eigs(op, n_eigvals, ncv);
- eigs.init();
- nconv = eigs.compute(1000, tol);
- eigval = eigs.eigenvalues();
- eigvec = eigs.eigenvectors();
- }
- else
- if(form_val == form_sm)
- {
- newarp::GenEigsSolver< T, newarp::EigsSelect::SMALLEST_MAGN, newarp::SparseGenMatProd<T> > eigs(op, n_eigvals, ncv);
- eigs.init();
- nconv = eigs.compute(1000, tol);
- eigval = eigs.eigenvalues();
- eigvec = eigs.eigenvectors();
- }
- else
- if(form_val == form_lr)
- {
- newarp::GenEigsSolver< T, newarp::EigsSelect::LARGEST_REAL, newarp::SparseGenMatProd<T> > eigs(op, n_eigvals, ncv);
- eigs.init();
- nconv = eigs.compute(1000, tol);
- eigval = eigs.eigenvalues();
- eigvec = eigs.eigenvectors();
- }
- else
- if(form_val == form_sr)
- {
- newarp::GenEigsSolver< T, newarp::EigsSelect::SMALLEST_REAL, newarp::SparseGenMatProd<T> > eigs(op, n_eigvals, ncv);
- eigs.init();
- nconv = eigs.compute(1000, tol);
- eigval = eigs.eigenvalues();
- eigvec = eigs.eigenvectors();
- }
- else
- if(form_val == form_li)
- {
- newarp::GenEigsSolver< T, newarp::EigsSelect::LARGEST_IMAG, newarp::SparseGenMatProd<T> > eigs(op, n_eigvals, ncv);
- eigs.init();
- nconv = eigs.compute(1000, tol);
- eigval = eigs.eigenvalues();
- eigvec = eigs.eigenvectors();
- }
- else
- if(form_val == form_si)
- {
- newarp::GenEigsSolver< T, newarp::EigsSelect::SMALLEST_IMAG, newarp::SparseGenMatProd<T> > eigs(op, n_eigvals, ncv);
- eigs.init();
- nconv = eigs.compute(1000, tol);
- eigval = eigs.eigenvalues();
- eigvec = eigs.eigenvectors();
- }
- }
- catch(const std::runtime_error&)
- {
- status = false;
- }
-
- if(status == true)
- {
- if(nconv == 0) { status = false; }
- }
-
- return status;
- }
- #else
- {
- arma_ignore(eigval);
- arma_ignore(eigvec);
- arma_ignore(X);
- arma_ignore(n_eigvals);
- arma_ignore(form_str);
- arma_ignore(default_tol);
-
- return false;
- }
- #endif
- }
- template<typename T>
- inline
- bool
- 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)
- {
- arma_extra_debug_sigprint();
-
- #if defined(ARMA_USE_ARPACK)
- {
- const form_type form_val = sp_auxlib::interpret_form_str(form_str);
-
- arma_debug_check( (form_val == form_none), "eigs_gen(): unknown form specified" );
-
- char which_lm[3] = "LM";
- char which_sm[3] = "SM";
- char which_lr[3] = "LR";
- char which_sr[3] = "SR";
- char which_li[3] = "LI";
- char which_si[3] = "SI";
-
- char* which;
-
- switch(form_val)
- {
- case form_lm: which = which_lm; break;
- case form_sm: which = which_sm; break;
- case form_lr: which = which_lr; break;
- case form_sr: which = which_sr; break;
- case form_li: which = which_li; break;
- case form_si: which = which_si; break;
-
- default: which = which_lm;
- }
-
-
- // Make sure it's square.
- arma_debug_check( (X.n_rows != X.n_cols), "eigs_gen(): given matrix must be square sized" );
-
- // Make sure we aren't asking for every eigenvalue.
- 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" );
-
- // If the matrix is empty, the case is trivial.
- if( (X.n_cols == 0) || (n_eigvals == 0) ) // We already know n_cols == n_rows.
- {
- eigval.reset();
- eigvec.reset();
- return true;
- }
-
- // Set up variables that get used for neupd().
- blas_int n, ncv, ldv, lworkl, info;
- T tol = default_tol;
- podarray<T> resid, v, workd, workl;
- podarray<blas_int> iparam, ipntr;
- podarray<T> rwork; // Not used in the real case.
-
- run_aupd(n_eigvals, which, X, false /* gen, not sym */, n, tol, resid, ncv, v, ldv, iparam, ipntr, workd, workl, lworkl, rwork, info);
-
- if(info != 0)
- {
- return false;
- }
- // The process has converged, and now we need to recover the actual eigenvectors using neupd().
- blas_int rvec = 1; // .TRUE
- blas_int nev = n_eigvals;
-
- char howmny = 'A';
- char bmat = 'I'; // We are considering the standard eigenvalue problem.
-
- podarray<blas_int> select(ncv); // Logical array of dimension NCV.
- podarray<T> dr(nev + 1); // Real array of dimension NEV + 1.
- podarray<T> di(nev + 1); // Real array of dimension NEV + 1.
- podarray<T> z(n * (nev + 1)); // Real N by NEV array if HOWMNY = 'A'.
- blas_int ldz = n;
- podarray<T> workev(3 * ncv);
-
- dr.zeros();
- di.zeros();
- z.zeros();
-
- 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);
-
- // Check for errors.
- if(info != 0)
- {
- arma_debug_warn("eigs_gen(): ARPACK error ", info, " in neupd()");
- return false;
- }
-
- // Put it into the outputs.
- eigval.set_size(n_eigvals);
- eigvec.zeros(n, n_eigvals);
-
- for (uword i = 0; i < n_eigvals; ++i)
- {
- eigval[i] = std::complex<T>(dr[i], di[i]);
- }
-
- // Now recover the eigenvectors.
- for (uword i = 0; i < n_eigvals; ++i)
- {
- // ARPACK ?neupd lays things out kinda odd in memory;
- // so does LAPACK ?geev -- see auxlib::eig_gen()
- if((i < n_eigvals - 1) && (eigval[i] == std::conj(eigval[i + 1])))
- {
- for (uword j = 0; j < uword(n); ++j)
- {
- eigvec.at(j, i) = std::complex<T>(z[n * i + j], z[n * (i + 1) + j]);
- eigvec.at(j, i + 1) = std::complex<T>(z[n * i + j], -z[n * (i + 1) + j]);
- }
- ++i; // Skip the next one.
- }
- else
- if((i == n_eigvals - 1) && (std::complex<T>(eigval[i]).imag() != 0.0))
- {
- // We don't have the matched conjugate eigenvalue.
- for (uword j = 0; j < uword(n); ++j)
- {
- eigvec.at(j, i) = std::complex<T>(z[n * i + j], z[n * (i + 1) + j]);
- }
- }
- else
- {
- // The eigenvector is entirely real.
- for (uword j = 0; j < uword(n); ++j)
- {
- eigvec.at(j, i) = std::complex<T>(z[n * i + j], T(0));
- }
- }
- }
-
- return (info == 0);
- }
- #else
- {
- arma_ignore(eigval);
- arma_ignore(eigvec);
- arma_ignore(X);
- arma_ignore(n_eigvals);
- arma_ignore(form_str);
- arma_ignore(default_tol);
-
- return false;
- }
- #endif
- }
- //! immediate eigendecomposition of non-symmetric complex sparse object
- template<typename T, typename T1>
- inline
- bool
- 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)
- {
- arma_extra_debug_sigprint();
-
- #if defined(ARMA_USE_ARPACK)
- {
- typedef typename std::complex<T> eT;
-
- const form_type form_val = sp_auxlib::interpret_form_str(form_str);
-
- arma_debug_check( (form_val == form_none), "eigs_gen(): unknown form specified" );
-
- char which_lm[3] = "LM";
- char which_sm[3] = "SM";
- char which_lr[3] = "LR";
- char which_sr[3] = "SR";
- char which_li[3] = "LI";
- char which_si[3] = "SI";
-
- char* which;
-
- switch(form_val)
- {
- case form_lm: which = which_lm; break;
- case form_sm: which = which_sm; break;
- case form_lr: which = which_lr; break;
- case form_sr: which = which_sr; break;
- case form_li: which = which_li; break;
- case form_si: which = which_si; break;
-
- default: which = which_lm;
- }
-
- const unwrap_spmat<T1> U(X_expr.get_ref());
-
- const SpMat<eT>& X = U.M;
-
- // Make sure it's square.
- arma_debug_check( (X.n_rows != X.n_cols), "eigs_gen(): given matrix must be square sized" );
-
- // Make sure we aren't asking for every eigenvalue.
- 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" );
-
- // If the matrix is empty, the case is trivial.
- if( (X.n_cols == 0) || (n_eigvals == 0) ) // We already know n_cols == n_rows.
- {
- eigval.reset();
- eigvec.reset();
- return true;
- }
-
- // Set up variables that get used for neupd().
- blas_int n, ncv, ldv, lworkl, info;
- T tol = default_tol;
- podarray< std::complex<T> > resid, v, workd, workl;
- podarray<blas_int> iparam, ipntr;
- podarray<T> rwork;
-
- run_aupd(n_eigvals, which, X, false /* gen, not sym */, n, tol, resid, ncv, v, ldv, iparam, ipntr, workd, workl, lworkl, rwork, info);
-
- if(info != 0)
- {
- return false;
- }
-
- // The process has converged, and now we need to recover the actual eigenvectors using neupd().
- blas_int rvec = 1; // .TRUE
- blas_int nev = n_eigvals;
-
- char howmny = 'A';
- char bmat = 'I'; // We are considering the standard eigenvalue problem.
-
- podarray<blas_int> select(ncv); // Logical array of dimension NCV.
- podarray<std::complex<T> > d(nev + 1); // Real array of dimension NEV + 1.
- podarray<std::complex<T> > z(n * nev); // Real N by NEV array if HOWMNY = 'A'.
- blas_int ldz = n;
- podarray<std::complex<T> > workev(2 * ncv);
-
- // Prepare the outputs; neupd() will write directly to them.
- eigval.zeros(n_eigvals);
- eigvec.zeros(n, n_eigvals);
- std::complex<T> sigma;
-
- arpack::neupd(&rvec, &howmny, select.memptr(), eigval.memptr(),
- (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);
-
- // Check for errors.
- if(info != 0)
- {
- arma_debug_warn("eigs_gen(): ARPACK error ", info, " in neupd()");
- return false;
- }
-
- return (info == 0);
- }
- #else
- {
- arma_ignore(eigval);
- arma_ignore(eigvec);
- arma_ignore(X_expr);
- arma_ignore(n_eigvals);
- arma_ignore(form_str);
- arma_ignore(default_tol);
-
- arma_stop_logic_error("eigs_gen(): use of ARPACK must be enabled for decomposition of complex matrices");
- return false;
- }
- #endif
- }
- template<typename T1, typename T2>
- inline
- bool
- 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)
- {
- arma_extra_debug_sigprint();
-
- #if defined(ARMA_USE_SUPERLU)
- {
- typedef typename T1::elem_type eT;
-
- superlu::superlu_options_t options;
- sp_auxlib::set_superlu_opts(options, user_opts);
-
- const unwrap_spmat<T1> tmp1(A_expr.get_ref());
- const SpMat<eT>& A = tmp1.M;
-
- X = B_expr.get_ref(); // superlu::gssv() uses X as input (the B matrix) and as output (the solution)
-
- if(A.n_rows > A.n_cols)
- {
- arma_stop_logic_error("spsolve(): solving over-determined systems currently not supported");
- X.soft_reset();
- return false;
- }
- else if(A.n_rows < A.n_cols)
- {
- arma_stop_logic_error("spsolve(): solving under-determined systems currently not supported");
- X.soft_reset();
- return false;
- }
-
- arma_debug_check( (A.n_rows != X.n_rows), "spsolve(): number of rows in the given objects must be the same" );
-
- if(A.is_empty() || X.is_empty())
- {
- X.zeros(A.n_cols, X.n_cols);
- return true;
- }
-
- if(A.n_nonzero == uword(0)) { X.soft_reset(); return false; }
-
- if(arma_config::debug)
- {
- bool overflow;
-
- overflow = (A.n_nonzero > INT_MAX);
- overflow = (A.n_rows > INT_MAX) || overflow;
- overflow = (A.n_cols > INT_MAX) || overflow;
- overflow = (X.n_rows > INT_MAX) || overflow;
- overflow = (X.n_cols > INT_MAX) || overflow;
-
- if(overflow)
- {
- arma_stop_runtime_error("spsolve(): integer overflow: matrix dimensions are too large for integer type used by SuperLU");
- return false;
- }
- }
-
- superlu::SuperMatrix x; arrayops::inplace_set(reinterpret_cast<char*>(&x), char(0), sizeof(superlu::SuperMatrix));
- superlu::SuperMatrix a; arrayops::inplace_set(reinterpret_cast<char*>(&a), char(0), sizeof(superlu::SuperMatrix));
-
- const bool status_x = wrap_to_supermatrix(x, X);
- const bool status_a = copy_to_supermatrix(a, A);
-
- if( (status_x == false) || (status_a == false) )
- {
- destroy_supermatrix(a);
- destroy_supermatrix(x);
- X.soft_reset();
- return false;
- }
-
- superlu::SuperMatrix l; arrayops::inplace_set(reinterpret_cast<char*>(&l), char(0), sizeof(superlu::SuperMatrix));
- superlu::SuperMatrix u; arrayops::inplace_set(reinterpret_cast<char*>(&u), char(0), sizeof(superlu::SuperMatrix));
-
- // paranoia: use SuperLU's memory allocation, in case it reallocs
-
- int* perm_c = (int*) superlu::malloc( (A.n_cols+1) * sizeof(int)); // extra paranoia: increase array length by 1
- int* perm_r = (int*) superlu::malloc( (A.n_rows+1) * sizeof(int));
-
- arma_check_bad_alloc( (perm_c == 0), "spsolve(): out of memory" );
- arma_check_bad_alloc( (perm_r == 0), "spsolve(): out of memory" );
-
- arrayops::inplace_set(perm_c, 0, A.n_cols+1);
- arrayops::inplace_set(perm_r, 0, A.n_rows+1);
-
- superlu::SuperLUStat_t stat;
- superlu::init_stat(&stat);
-
- int info = 0; // Return code.
-
- arma_extra_debug_print("superlu::gssv()");
- superlu::gssv<eT>(&options, &a, perm_c, perm_r, &l, &u, &x, &stat, &info);
-
-
- // Process the return code.
- if( (info > 0) && (info <= int(A.n_cols)) )
- {
- // std::ostringstream tmp;
- // tmp << "spsolve(): could not solve system; LU factorisation completed, but detected zero in U(" << (info-1) << ',' << (info-1) << ')';
- // arma_debug_warn(tmp.str());
- }
- else
- if(info > int(A.n_cols))
- {
- arma_debug_warn("spsolve(): memory allocation failure: could not allocate ", (info - int(A.n_cols)), " bytes");
- }
- else
- if(info < 0)
- {
- arma_debug_warn("spsolve(): unknown SuperLU error code from gssv(): ", info);
- }
-
-
- superlu::free_stat(&stat);
-
- superlu::free(perm_c);
- superlu::free(perm_r);
-
- destroy_supermatrix(u);
- destroy_supermatrix(l);
- destroy_supermatrix(a);
- destroy_supermatrix(x); // No need to extract the data from x, since it's using the same memory as X
-
- return (info == 0);
- }
- #else
- {
- arma_ignore(X);
- arma_ignore(A_expr);
- arma_ignore(B_expr);
- arma_ignore(user_opts);
- arma_stop_logic_error("spsolve(): use of SuperLU must be enabled");
- return false;
- }
- #endif
- }
- template<typename T1, typename T2>
- inline
- bool
- 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)
- {
- arma_extra_debug_sigprint();
-
- #if defined(ARMA_USE_SUPERLU)
- {
- typedef typename T1::pod_type T;
- typedef typename T1::elem_type eT;
-
- superlu::superlu_options_t options;
- sp_auxlib::set_superlu_opts(options, user_opts);
-
- const unwrap_spmat<T1> tmp1(A_expr.get_ref());
- const SpMat<eT>& A = tmp1.M;
-
- const unwrap<T2> tmp2(B_expr.get_ref());
- const Mat<eT>& B_unwrap = tmp2.M;
-
- const bool B_is_modified = ( (user_opts.equilibrate) || (&B_unwrap == &X) );
-
- Mat<eT> B_copy; if(B_is_modified) { B_copy = B_unwrap; }
-
- const Mat<eT>& B = (B_is_modified) ? B_copy : B_unwrap;
-
- if(A.n_rows > A.n_cols)
- {
- arma_stop_logic_error("spsolve(): solving over-determined systems currently not supported");
- X.soft_reset();
- return false;
- }
- else if(A.n_rows < A.n_cols)
- {
- arma_stop_logic_error("spsolve(): solving under-determined systems currently not supported");
- X.soft_reset();
- return false;
- }
-
- arma_debug_check( (A.n_rows != B.n_rows), "spsolve(): number of rows in the given objects must be the same" );
-
- X.zeros(A.n_cols, B.n_cols); // set the elements to zero, as we don't trust the SuperLU spaghetti code
-
- if(A.is_empty() || B.is_empty())
- {
- return true;
- }
-
- if(A.n_nonzero == uword(0)) { X.soft_reset(); return false; }
-
- if(arma_config::debug)
- {
- bool overflow;
-
- overflow = (A.n_nonzero > INT_MAX);
- overflow = (A.n_rows > INT_MAX) || overflow;
- overflow = (A.n_cols > INT_MAX) || overflow;
- overflow = (B.n_rows > INT_MAX) || overflow;
- overflow = (B.n_cols > INT_MAX) || overflow;
- overflow = (X.n_rows > INT_MAX) || overflow;
- overflow = (X.n_cols > INT_MAX) || overflow;
-
- if(overflow)
- {
- arma_stop_runtime_error("spsolve(): integer overflow: matrix dimensions are too large for integer type used by SuperLU");
- return false;
- }
- }
-
- superlu::SuperMatrix x; arrayops::inplace_set(reinterpret_cast<char*>(&x), char(0), sizeof(superlu::SuperMatrix));
- superlu::SuperMatrix a; arrayops::inplace_set(reinterpret_cast<char*>(&a), char(0), sizeof(superlu::SuperMatrix));
- superlu::SuperMatrix b; arrayops::inplace_set(reinterpret_cast<char*>(&b), char(0), sizeof(superlu::SuperMatrix));
-
- const bool status_x = wrap_to_supermatrix(x, X);
- const bool status_a = copy_to_supermatrix(a, A); // NOTE: superlu::gssvx() modifies 'a' if equilibration is enabled
- const bool status_b = wrap_to_supermatrix(b, B); // NOTE: superlu::gssvx() modifies 'b' if equilibration is enabled
-
- if( (status_x == false) || (status_a == false) || (status_b == false) )
- {
- destroy_supermatrix(x);
- destroy_supermatrix(a);
- destroy_supermatrix(b);
- X.soft_reset();
- return false;
- }
-
- superlu::SuperMatrix l; arrayops::inplace_set(reinterpret_cast<char*>(&l), char(0), sizeof(superlu::SuperMatrix));
- superlu::SuperMatrix u; arrayops::inplace_set(reinterpret_cast<char*>(&u), char(0), sizeof(superlu::SuperMatrix));
-
- // paranoia: use SuperLU's memory allocation, in case it reallocs
-
- int* perm_c = (int*) superlu::malloc( (A.n_cols+1) * sizeof(int) ); // extra paranoia: increase array length by 1
- int* perm_r = (int*) superlu::malloc( (A.n_rows+1) * sizeof(int) );
- int* etree = (int*) superlu::malloc( (A.n_cols+1) * sizeof(int) );
-
- T* R = (T*) superlu::malloc( (A.n_rows+1) * sizeof(T) );
- T* C = (T*) superlu::malloc( (A.n_cols+1) * sizeof(T) );
- T* ferr = (T*) superlu::malloc( (B.n_cols+1) * sizeof(T) );
- T* berr = (T*) superlu::malloc( (B.n_cols+1) * sizeof(T) );
-
- arma_check_bad_alloc( (perm_c == 0), "spsolve(): out of memory" );
- arma_check_bad_alloc( (perm_r == 0), "spsolve(): out of memory" );
- arma_check_bad_alloc( (etree == 0), "spsolve(): out of memory" );
-
- arma_check_bad_alloc( (R == 0), "spsolve(): out of memory" );
- arma_check_bad_alloc( (C == 0), "spsolve(): out of memory" );
- arma_check_bad_alloc( (ferr == 0), "spsolve(): out of memory" );
- arma_check_bad_alloc( (berr == 0), "spsolve(): out of memory" );
-
- arrayops::inplace_set(perm_c, int(0), A.n_cols+1);
- arrayops::inplace_set(perm_r, int(0), A.n_rows+1);
- arrayops::inplace_set(etree, int(0), A.n_cols+1);
-
- arrayops::inplace_set(R, T(0), A.n_rows+1);
- arrayops::inplace_set(C, T(0), A.n_cols+1);
- arrayops::inplace_set(ferr, T(0), B.n_cols+1);
- arrayops::inplace_set(berr, T(0), B.n_cols+1);
-
- superlu::GlobalLU_t glu;
- arrayops::inplace_set(reinterpret_cast<char*>(&glu), char(0), sizeof(superlu::GlobalLU_t));
-
- superlu::mem_usage_t mu;
- arrayops::inplace_set(reinterpret_cast<char*>(&mu), char(0), sizeof(superlu::mem_usage_t));
-
- superlu::SuperLUStat_t stat;
- superlu::init_stat(&stat);
-
- char equed[8]; // extra characters for paranoia
- T rpg = T(0);
- T rcond = T(0);
- int info = int(0); // Return code.
-
- char work[8];
- int lwork = int(0); // 0 means superlu will allocate memory
-
- arma_extra_debug_print("superlu::gssvx()");
- 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);
-
- bool status = false;
-
- // Process the return code.
- if(info == 0)
- {
- status = true;
- }
- if( (info > 0) && (info <= int(A.n_cols)) )
- {
- // std::ostringstream tmp;
- // tmp << "spsolve(): could not solve system; LU factorisation completed, but detected zero in U(" << (info-1) << ',' << (info-1) << ')';
- // arma_debug_warn(tmp.str());
- }
- else
- if( (info == int(A.n_cols+1)) && (user_opts.allow_ugly) )
- {
- arma_debug_warn("spsolve(): system is singular to working precision (rcond: ", rcond, ")");
- status = true;
- }
- else
- if(info > int(A.n_cols+1))
- {
- arma_debug_warn("spsolve(): memory allocation failure: could not allocate ", (info - int(A.n_cols)), " bytes");
- }
- else
- if(info < 0)
- {
- arma_debug_warn("spsolve(): unknown SuperLU error code from gssvx(): ", info);
- }
-
- superlu::free_stat(&stat);
-
- superlu::free(berr);
- superlu::free(ferr);
- superlu::free(C);
- superlu::free(R);
- superlu::free(etree);
- superlu::free(perm_r);
- superlu::free(perm_c);
-
- destroy_supermatrix(u);
- destroy_supermatrix(l);
- destroy_supermatrix(b);
- destroy_supermatrix(a);
- destroy_supermatrix(x); // No need to extract the data from x, since it's using the same memory as X
-
- out_rcond = rcond;
-
- return status;
- }
- #else
- {
- arma_ignore(X);
- arma_ignore(out_rcond);
- arma_ignore(A_expr);
- arma_ignore(B_expr);
- arma_ignore(user_opts);
- arma_stop_logic_error("spsolve(): use of SuperLU must be enabled");
- return false;
- }
- #endif
- }
- #if defined(ARMA_USE_SUPERLU)
-
- inline
- void
- sp_auxlib::set_superlu_opts(superlu::superlu_options_t& options, const superlu_opts& user_opts)
- {
- arma_extra_debug_sigprint();
-
- // default options as the starting point
- superlu::set_default_opts(&options);
-
- // our settings
- options.Trans = superlu::NOTRANS;
- options.ConditionNumber = superlu::YES;
-
- // process user_opts
-
- if(user_opts.equilibrate == true) { options.Equil = superlu::YES; }
- if(user_opts.equilibrate == false) { options.Equil = superlu::NO; }
-
- if(user_opts.symmetric == true) { options.SymmetricMode = superlu::YES; }
- if(user_opts.symmetric == false) { options.SymmetricMode = superlu::NO; }
-
- options.DiagPivotThresh = user_opts.pivot_thresh;
-
- if(user_opts.permutation == superlu_opts::NATURAL) { options.ColPerm = superlu::NATURAL; }
- if(user_opts.permutation == superlu_opts::MMD_ATA) { options.ColPerm = superlu::MMD_ATA; }
- if(user_opts.permutation == superlu_opts::MMD_AT_PLUS_A) { options.ColPerm = superlu::MMD_AT_PLUS_A; }
- if(user_opts.permutation == superlu_opts::COLAMD) { options.ColPerm = superlu::COLAMD; }
-
- if(user_opts.refine == superlu_opts::REF_NONE) { options.IterRefine = superlu::NOREFINE; }
- if(user_opts.refine == superlu_opts::REF_SINGLE) { options.IterRefine = superlu::SLU_SINGLE; }
- if(user_opts.refine == superlu_opts::REF_DOUBLE) { options.IterRefine = superlu::SLU_DOUBLE; }
- if(user_opts.refine == superlu_opts::REF_EXTRA) { options.IterRefine = superlu::SLU_EXTRA; }
- }
-
-
-
- template<typename eT>
- inline
- bool
- sp_auxlib::copy_to_supermatrix(superlu::SuperMatrix& out, const SpMat<eT>& A)
- {
- arma_extra_debug_sigprint();
-
- // We store in column-major CSC.
- out.Stype = superlu::SLU_NC;
-
- if(is_float<eT>::value)
- {
- out.Dtype = superlu::SLU_S;
- }
- else
- if(is_double<eT>::value)
- {
- out.Dtype = superlu::SLU_D;
- }
- else
- if(is_cx_float<eT>::value)
- {
- out.Dtype = superlu::SLU_C;
- }
- else
- if(is_cx_double<eT>::value)
- {
- out.Dtype = superlu::SLU_Z;
- }
-
- out.Mtype = superlu::SLU_GE; // Just a general matrix. We don't know more now.
-
- // We have to actually create the object which stores the data.
- // This gets cleaned by destroy_supermatrix().
- // We have to use SuperLU's problematic memory allocation routines since they are
- // not guaranteed to be new and delete. See the comments in def_superlu.hpp
- superlu::NCformat* nc = (superlu::NCformat*)superlu::malloc(sizeof(superlu::NCformat));
-
- if(nc == NULL) { return false; }
-
- A.sync();
-
- nc->nnz = A.n_nonzero;
- nc->nzval = (void*) superlu::malloc(sizeof(eT) * A.n_nonzero );
- nc->colptr = (superlu::int_t*)superlu::malloc(sizeof(superlu::int_t) * (A.n_cols + 1));
- nc->rowind = (superlu::int_t*)superlu::malloc(sizeof(superlu::int_t) * A.n_nonzero );
-
- if( (nc->nzval == NULL) || (nc->colptr == NULL) || (nc->rowind == NULL) ) { return false; }
-
- // Fill the matrix.
- arrayops::copy((eT*) nc->nzval, A.values, A.n_nonzero);
-
- // // These have to be copied by hand, because the types may differ.
- // for (uword i = 0; i <= A.n_cols; ++i) { nc->colptr[i] = (int_t) A.col_ptrs[i]; }
- // for (uword i = 0; i < A.n_nonzero; ++i) { nc->rowind[i] = (int_t) A.row_indices[i]; }
-
- arrayops::convert(nc->colptr, A.col_ptrs, A.n_cols+1 );
- arrayops::convert(nc->rowind, A.row_indices, A.n_nonzero);
-
- out.nrow = A.n_rows;
- out.ncol = A.n_cols;
- out.Store = (void*) nc;
-
- return true;
- }
-
-
-
- template<typename eT>
- inline
- bool
- sp_auxlib::wrap_to_supermatrix(superlu::SuperMatrix& out, const Mat<eT>& A)
- {
- arma_extra_debug_sigprint();
-
- // NOTE: this function re-uses memory from matrix A
-
- // This is being stored as a dense matrix.
- out.Stype = superlu::SLU_DN;
-
- if(is_float<eT>::value)
- {
- out.Dtype = superlu::SLU_S;
- }
- else
- if(is_double<eT>::value)
- {
- out.Dtype = superlu::SLU_D;
- }
- else
- if(is_cx_float<eT>::value)
- {
- out.Dtype = superlu::SLU_C;
- }
- else
- if(is_cx_double<eT>::value)
- {
- out.Dtype = superlu::SLU_Z;
- }
-
- out.Mtype = superlu::SLU_GE;
-
- // We have to create the object that stores the data.
- superlu::DNformat* dn = (superlu::DNformat*)superlu::malloc(sizeof(superlu::DNformat));
-
- if(dn == NULL) { return false; }
-
- dn->lda = A.n_rows;
- dn->nzval = (void*) A.memptr(); // re-use memory instead of copying
-
- out.nrow = A.n_rows;
- out.ncol = A.n_cols;
- out.Store = (void*) dn;
-
- return true;
- }
-
-
-
- inline
- void
- sp_auxlib::destroy_supermatrix(superlu::SuperMatrix& out)
- {
- arma_extra_debug_sigprint();
-
- // Clean up.
- if(out.Stype == superlu::SLU_NC)
- {
- superlu::destroy_compcol_mat(&out);
- }
- else
- if(out.Stype == superlu::SLU_DN)
- {
- // superlu::destroy_dense_mat(&out);
-
- // since dn->nzval is set to re-use memory from a Mat object (which manages its own memory),
- // we cannot simply call superlu::destroy_dense_mat().
- // Only the out.Store structure can be freed.
-
- superlu::DNformat* dn = (superlu::DNformat*) out.Store;
-
- if(dn != NULL) { superlu::free(dn); }
- }
- else
- if(out.Stype == superlu::SLU_SC)
- {
- superlu::destroy_supernode_mat(&out);
- }
- else
- {
- // Uh, crap.
-
- std::ostringstream tmp;
-
- tmp << "sp_auxlib::destroy_supermatrix(): unhandled Stype" << std::endl;
- tmp << "Stype val: " << out.Stype << std::endl;
- tmp << "Stype name: ";
-
- if(out.Stype == superlu::SLU_NC) { tmp << "SLU_NC"; }
- if(out.Stype == superlu::SLU_NCP) { tmp << "SLU_NCP"; }
- if(out.Stype == superlu::SLU_NR) { tmp << "SLU_NR"; }
- if(out.Stype == superlu::SLU_SC) { tmp << "SLU_SC"; }
- if(out.Stype == superlu::SLU_SCP) { tmp << "SLU_SCP"; }
- if(out.Stype == superlu::SLU_SR) { tmp << "SLU_SR"; }
- if(out.Stype == superlu::SLU_DN) { tmp << "SLU_DN"; }
- if(out.Stype == superlu::SLU_NR_loc) { tmp << "SLU_NR_loc"; }
-
- arma_debug_warn(tmp.str());
- arma_stop_runtime_error("internal error: sp_auxlib::destroy_supermatrix()");
- }
- }
-
- #endif
- template<typename eT, typename T>
- inline
- void
- sp_auxlib::run_aupd
- (
- const uword n_eigvals, char* which, const SpMat<T>& X, const bool sym,
- blas_int& n, eT& tol,
- podarray<T>& resid, blas_int& ncv, podarray<T>& v, blas_int& ldv,
- podarray<blas_int>& iparam, podarray<blas_int>& ipntr,
- podarray<T>& workd, podarray<T>& workl, blas_int& lworkl, podarray<eT>& rwork,
- blas_int& info
- )
- {
- #if defined(ARMA_USE_ARPACK)
- {
- // ARPACK provides a "reverse communication interface" which is an
- // entertainingly archaic FORTRAN software engineering technique that
- // basically means that we call saupd()/naupd() and it tells us with some
- // return code what we need to do next (usually a matrix-vector product) and
- // then call it again. So this results in some type of iterative process
- // where we call saupd()/naupd() many times.
- blas_int ido = 0; // This must be 0 for the first call.
- char bmat = 'I'; // We are considering the standard eigenvalue problem.
- n = X.n_rows; // The size of the matrix.
- blas_int nev = n_eigvals;
-
- resid.set_size(n);
-
- // Two contraints on NCV: (NCV > NEV + 2) and (NCV <= N)
- //
- // We're calling either arpack::saupd() or arpack::naupd(),
- // which have slighly different minimum constraint and recommended value for NCV:
- // http://www.caam.rice.edu/software/ARPACK/UG/node136.html
- // http://www.caam.rice.edu/software/ARPACK/UG/node138.html
-
- ncv = nev + 2 + 1;
-
- if (ncv < (2 * nev + 1)) { ncv = 2 * nev + 1; }
- if (ncv > n ) { ncv = n; }
-
- v.set_size(n * ncv); // Array N by NCV (output).
- rwork.set_size(ncv); // Work array of size NCV for complex calls.
- ldv = n; // "Leading dimension of V exactly as declared in the calling program."
-
- // IPARAM: integer array of length 11.
- iparam.zeros(11);
- iparam(0) = 1; // Exact shifts (not provided by us).
- iparam(2) = 1000; // Maximum iterations; all the examples use 300, but they were written in the ancient times.
- iparam(6) = 1; // Mode 1: A * x = lambda * x.
-
- // IPNTR: integer array of length 14 (output).
- ipntr.set_size(14);
-
- // Real work array used in the basic Arnoldi iteration for reverse communication.
- workd.set_size(3 * n);
-
- // lworkl must be at least 3 * NCV^2 + 6 * NCV.
- lworkl = 3 * (ncv * ncv) + 6 * ncv;
-
- // Real work array of length lworkl.
- workl.set_size(lworkl);
-
- info = 0; // Set to 0 initially to use random initial vector.
-
- // All the parameters have been set or created. Time to loop a lot.
- while (ido != 99)
- {
- // Call saupd() or naupd() with the current parameters.
- if(sym)
- arpack::saupd(&ido, &bmat, &n, which, &nev, &tol, resid.memptr(), &ncv, v.memptr(), &ldv, iparam.memptr(), ipntr.memptr(), workd.memptr(), workl.memptr(), &lworkl, &info);
- else
- 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);
-
- // What do we do now?
- switch (ido)
- {
- case -1:
- // fallthrough
- case 1:
- {
- // We need to calculate the matrix-vector multiplication y = OP * x
- // where x is of length n and starts at workd(ipntr(0)), and y is of
- // length n and starts at workd(ipntr(1)).
-
- // operator*(sp_mat, vec) doesn't properly put the result into the
- // right place so we'll just reimplement it here for now...
-
- // Set the output to point at the right memory. We have to subtract
- // one from FORTRAN pointers...
- Col<T> out(workd.memptr() + ipntr(1) - 1, n, false /* don't copy */);
- // Set the input to point at the right memory.
- Col<T> in(workd.memptr() + ipntr(0) - 1, n, false /* don't copy */);
-
- out.zeros();
- typename SpMat<T>::const_iterator x_it = X.begin();
- typename SpMat<T>::const_iterator x_it_end = X.end();
-
- while(x_it != x_it_end)
- {
- out[x_it.row()] += (*x_it) * in[x_it.col()];
- ++x_it;
- }
-
- // No need to modify memory further since it was all done in-place.
-
- break;
- }
- case 99:
- // Nothing to do here, things have converged.
- break;
- default:
- {
- return; // Parent frame can look at the value of info.
- }
- }
- }
-
- // The process has ended; check the return code.
- if( (info != 0) && (info != 1) )
- {
- // Print warnings if there was a failure.
-
- if(sym)
- {
- arma_debug_warn("eigs_sym(): ARPACK error ", info, " in saupd()");
- }
- else
- {
- arma_debug_warn("eigs_gen(): ARPACK error ", info, " in naupd()");
- }
-
- return; // Parent frame can look at the value of info.
- }
- }
- #else
- arma_ignore(n_eigvals);
- arma_ignore(which);
- arma_ignore(X);
- arma_ignore(sym);
- arma_ignore(n);
- arma_ignore(tol);
- arma_ignore(resid);
- arma_ignore(ncv);
- arma_ignore(v);
- arma_ignore(ldv);
- arma_ignore(iparam);
- arma_ignore(ipntr);
- arma_ignore(workd);
- arma_ignore(workl);
- arma_ignore(lworkl);
- arma_ignore(rwork);
- arma_ignore(info);
- #endif
- }
- template<typename eT>
- inline
- bool
- sp_auxlib::rudimentary_sym_check(const SpMat<eT>& X)
- {
- arma_extra_debug_sigprint();
-
- if(X.n_rows != X.n_cols) { return false; }
-
- const eT tol = eT(10000) * std::numeric_limits<eT>::epsilon(); // allow some leeway
-
- typename SpMat<eT>::const_iterator it = X.begin();
- typename SpMat<eT>::const_iterator it_end = X.end();
-
- const uword n_check_limit = (std::max)( uword(2), uword(X.n_nonzero/100) );
-
- uword n_check = 1;
-
- while( (it != it_end) && (n_check <= n_check_limit) )
- {
- const uword it_row = it.row();
- const uword it_col = it.col();
-
- if(it_row != it_col)
- {
- const eT A = (*it);
- const eT B = X.at( it_col, it_row ); // deliberately swapped
-
- const eT C = (std::max)(std::abs(A), std::abs(B));
-
- const eT delta = std::abs(A - B);
-
- if(( (delta <= tol) || (delta <= (C * tol)) ) == false) { return false; }
-
- ++n_check;
- }
-
- ++it;
- }
-
- return true;
- }
- template<typename T>
- inline
- bool
- sp_auxlib::rudimentary_sym_check(const SpMat< std::complex<T> >& X)
- {
- arma_extra_debug_sigprint();
-
- // NOTE: the function name is a misnomer, as it checks for hermitian complex matrices;
- // NOTE: for simplicity of use, the function name is the same as for real matrices
-
- typedef typename std::complex<T> eT;
-
- if(X.n_rows != X.n_cols) { return false; }
-
- const T tol = T(10000) * std::numeric_limits<T>::epsilon(); // allow some leeway
-
- typename SpMat<eT>::const_iterator it = X.begin();
- typename SpMat<eT>::const_iterator it_end = X.end();
-
- const uword n_check_limit = (std::max)( uword(2), uword(X.n_nonzero/100) );
-
- uword n_check = 1;
-
- while( (it != it_end) && (n_check <= n_check_limit) )
- {
- const uword it_row = it.row();
- const uword it_col = it.col();
-
- if(it_row != it_col)
- {
- const eT A = (*it);
- const eT B = X.at( it_col, it_row ); // deliberately swapped
-
- const T C_real = (std::max)(std::abs(A.real()), std::abs(B.real()));
- const T C_imag = (std::max)(std::abs(A.imag()), std::abs(B.imag()));
-
- const T delta_real = std::abs(A.real() - B.real());
- const T delta_imag = std::abs(A.imag() + B.imag()); // take into account the conjugate
-
- const bool okay_real = ( (delta_real <= tol) || (delta_real <= (C_real * tol)) );
- const bool okay_imag = ( (delta_imag <= tol) || (delta_imag <= (C_imag * tol)) );
-
- if( (okay_real == false) || (okay_imag == false) ) { return false; }
-
- ++n_check;
- }
- else
- {
- const eT A = (*it);
-
- if(std::abs(A.imag()) > tol) { return false; }
- }
-
- ++it;
- }
-
- return true;
- }
- //! @}
|