eigs.cc

Go to the documentation of this file.
00001 /*
00002 
00003 Copyright (C) 2005-2012 David Bateman
00004 
00005 This file is part of Octave.
00006 
00007 Octave is free software; you can redistribute it and/or modify it
00008 under the terms of the GNU General Public License as published by the
00009 Free Software Foundation; either version 3 of the License, or (at your
00010 option) any later version.
00011 
00012 Octave is distributed in the hope that it will be useful, but WITHOUT
00013 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00014 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00015 for more details.
00016 
00017 You should have received a copy of the GNU General Public License
00018 along with Octave; see the file COPYING.  If not, see
00019 <http://www.gnu.org/licenses/>.
00020 
00021 */
00022 
00023 #ifdef HAVE_CONFIG_H
00024 #include <config.h>
00025 #endif
00026 
00027 #include "ov.h"
00028 #include "defun-dld.h"
00029 #include "error.h"
00030 #include "gripes.h"
00031 #include "quit.h"
00032 #include "variables.h"
00033 #include "ov-re-sparse.h"
00034 #include "ov-cx-sparse.h"
00035 #include "oct-map.h"
00036 #include "pager.h"
00037 #include "unwind-prot.h"
00038 
00039 #include "eigs-base.cc"
00040 
00041 // Global pointer for user defined function.
00042 static octave_function *eigs_fcn = 0;
00043 
00044 // Have we warned about imaginary values returned from user function?
00045 static bool warned_imaginary = false;
00046 
00047 // Is this a recursive call?
00048 static int call_depth = 0;
00049 
00050 ColumnVector
00051 eigs_func (const ColumnVector &x, int &eigs_error)
00052 {
00053   ColumnVector retval;
00054   octave_value_list args;
00055   args(0) = x;
00056 
00057   if (eigs_fcn)
00058     {
00059       octave_value_list tmp = eigs_fcn->do_multi_index_op (1, args);
00060 
00061       if (error_state)
00062         {
00063           eigs_error = 1;
00064           gripe_user_supplied_eval ("eigs");
00065           return retval;
00066         }
00067 
00068       if (tmp.length () && tmp(0).is_defined ())
00069         {
00070           if (! warned_imaginary && tmp(0).is_complex_type ())
00071             {
00072               warning ("eigs: ignoring imaginary part returned from user-supplied function");
00073               warned_imaginary = true;
00074             }
00075 
00076           retval = ColumnVector (tmp(0).vector_value ());
00077 
00078           if (error_state)
00079             {
00080               eigs_error = 1;
00081               gripe_user_supplied_eval ("eigs");
00082             }
00083         }
00084       else
00085         {
00086           eigs_error = 1;
00087           gripe_user_supplied_eval ("eigs");
00088         }
00089     }
00090 
00091   return retval;
00092 }
00093 
00094 ComplexColumnVector
00095 eigs_complex_func (const ComplexColumnVector &x, int &eigs_error)
00096 {
00097   ComplexColumnVector retval;
00098   octave_value_list args;
00099   args(0) = x;
00100 
00101   if (eigs_fcn)
00102     {
00103       octave_value_list tmp = eigs_fcn->do_multi_index_op (1, args);
00104 
00105       if (error_state)
00106         {
00107           eigs_error = 1;
00108           gripe_user_supplied_eval ("eigs");
00109           return retval;
00110         }
00111 
00112       if (tmp.length () && tmp(0).is_defined ())
00113         {
00114           retval = ComplexColumnVector (tmp(0).complex_vector_value ());
00115 
00116           if (error_state)
00117             {
00118               eigs_error = 1;
00119               gripe_user_supplied_eval ("eigs");
00120             }
00121         }
00122       else
00123         {
00124           eigs_error = 1;
00125           gripe_user_supplied_eval ("eigs");
00126         }
00127     }
00128 
00129   return retval;
00130 }
00131 
00132 DEFUN_DLD (eigs, args, nargout,
00133   "-*- texinfo -*-\n\
00134 @deftypefn  {Loadable Function} {@var{d} =} eigs (@var{A})\n\
00135 @deftypefnx {Loadable Function} {@var{d} =} eigs (@var{A}, @var{k})\n\
00136 @deftypefnx {Loadable Function} {@var{d} =} eigs (@var{A}, @var{k}, @var{sigma})\n\
00137 @deftypefnx {Loadable Function} {@var{d} =} eigs (@var{A}, @var{k}, @var{sigma}, @var{opts})\n\
00138 @deftypefnx {Loadable Function} {@var{d} =} eigs (@var{A}, @var{B})\n\
00139 @deftypefnx {Loadable Function} {@var{d} =} eigs (@var{A}, @var{B}, @var{k})\n\
00140 @deftypefnx {Loadable Function} {@var{d} =} eigs (@var{A}, @var{B}, @var{k}, @var{sigma})\n\
00141 @deftypefnx {Loadable Function} {@var{d} =} eigs (@var{A}, @var{B}, @var{k}, @var{sigma}, @var{opts})\n\
00142 @deftypefnx {Loadable Function} {@var{d} =} eigs (@var{af}, @var{n})\n\
00143 @deftypefnx {Loadable Function} {@var{d} =} eigs (@var{af}, @var{n}, @var{B})\n\
00144 @deftypefnx {Loadable Function} {@var{d} =} eigs (@var{af}, @var{n}, @var{k})\n\
00145 @deftypefnx {Loadable Function} {@var{d} =} eigs (@var{af}, @var{n}, @var{B}, @var{k})\n\
00146 @deftypefnx {Loadable Function} {@var{d} =} eigs (@var{af}, @var{n}, @var{k}, @var{sigma})\n\
00147 @deftypefnx {Loadable Function} {@var{d} =} eigs (@var{af}, @var{n}, @var{B}, @var{k}, @var{sigma})\n\
00148 @deftypefnx {Loadable Function} {@var{d} =} eigs (@var{af}, @var{n}, @var{k}, @var{sigma}, @var{opts})\n\
00149 @deftypefnx {Loadable Function} {@var{d} =} eigs (@var{af}, @var{n}, @var{B}, @var{k}, @var{sigma}, @var{opts})\n\
00150 @deftypefnx {Loadable Function} {[@var{V}, @var{d}] =} eigs (@var{A}, @dots{})\n\
00151 @deftypefnx {Loadable Function} {[@var{V}, @var{d}] =} eigs (@var{af}, @var{n}, @dots{})\n\
00152 @deftypefnx {Loadable Function} {[@var{V}, @var{d}, @var{flag}] =} eigs (@var{A}, @dots{})\n\
00153 @deftypefnx {Loadable Function} {[@var{V}, @var{d}, @var{flag}] =} eigs (@var{af}, @var{n}, @dots{})\n\
00154 Calculate a limited number of eigenvalues and eigenvectors of @var{A},\n\
00155 based on a selection criteria.  The number of eigenvalues and eigenvectors to\n\
00156 calculate is given by @var{k} and defaults to 6.\n\
00157 \n\
00158 By default, @code{eigs} solve the equation\n\
00159 @tex\n\
00160 $A \\nu = \\lambda \\nu$,\n\
00161 @end tex\n\
00162 @ifinfo\n\
00163 @code{A * v = lambda * v},\n\
00164 @end ifinfo\n\
00165 where\n\
00166 @tex\n\
00167 $\\lambda$ is a scalar representing one of the eigenvalues, and $\\nu$\n\
00168 @end tex\n\
00169 @ifinfo\n\
00170 @code{lambda} is a scalar representing one of the eigenvalues, and @code{v}\n\
00171 @end ifinfo\n\
00172 is the corresponding eigenvector.  If given the positive definite matrix\n\
00173 @var{B} then @code{eigs} solves the general eigenvalue equation\n\
00174 @tex\n\
00175 $A \\nu = \\lambda B \\nu$.\n\
00176 @end tex\n\
00177 @ifinfo\n\
00178 @code{A * v = lambda * B * v}.\n\
00179 @end ifinfo\n\
00180 \n\
00181 The argument @var{sigma} determines which eigenvalues are returned.\n\
00182 @var{sigma} can be either a scalar or a string.  When @var{sigma} is a\n\
00183 scalar, the @var{k} eigenvalues closest to @var{sigma} are returned.  If\n\
00184 @var{sigma} is a string, it must have one of the following values.\n\
00185 \n\
00186 @table @asis\n\
00187 @item 'lm'\n\
00188 Largest Magnitude (default).\n\
00189 \n\
00190 @item 'sm'\n\
00191 Smallest Magnitude.\n\
00192 \n\
00193 @item 'la'\n\
00194 Largest Algebraic (valid only for real symmetric problems).\n\
00195 \n\
00196 @item 'sa'\n\
00197 Smallest Algebraic (valid only for real symmetric problems).\n\
00198 \n\
00199 @item 'be'\n\
00200 Both Ends, with one more from the high-end if @var{k} is odd (valid only for\n\
00201 real symmetric problems).\n\
00202 \n\
00203 @item 'lr'\n\
00204 Largest Real part (valid only for complex or unsymmetric problems).\n\
00205 \n\
00206 @item 'sr'\n\
00207 Smallest Real part (valid only for complex or unsymmetric problems).\n\
00208 \n\
00209 @item 'li'\n\
00210 Largest Imaginary part (valid only for complex or unsymmetric problems).\n\
00211 \n\
00212 @item 'si'\n\
00213 Smallest Imaginary part (valid only for complex or unsymmetric problems).\n\
00214 @end table\n\
00215 \n\
00216 If @var{opts} is given, it is a structure defining possible options that\n\
00217 @code{eigs} should use.  The fields of the @var{opts} structure are:\n\
00218 \n\
00219 @table @code\n\
00220 @item issym\n\
00221 If @var{af} is given, then flags whether the function @var{af} defines a\n\
00222 symmetric problem.  It is ignored if @var{A} is given.  The default is false.\n\
00223 \n\
00224 @item isreal\n\
00225 If @var{af} is given, then flags whether the function @var{af} defines a\n\
00226 real problem.  It is ignored if @var{A} is given.  The default is true.\n\
00227 \n\
00228 @item tol\n\
00229 Defines the required convergence tolerance, calculated as\n\
00230 @code{tol * norm (A)}.  The default is @code{eps}.\n\
00231 \n\
00232 @item maxit\n\
00233 The maximum number of iterations.  The default is 300.\n\
00234 \n\
00235 @item p\n\
00236 The number of Lanzcos basis vectors to use.  More vectors will result in\n\
00237 faster convergence, but a greater use of memory.  The optimal value of\n\
00238 @code{p} is problem dependent and should be in the range @var{k} to @var{n}.\n\
00239 The default value is @code{2 * @var{k}}.\n\
00240 \n\
00241 @item v0\n\
00242 The starting vector for the algorithm.  An initial vector close to the\n\
00243 final vector will speed up convergence.  The default is for @sc{arpack}\n\
00244 to randomly generate a starting vector.  If specified, @code{v0} must be\n\
00245 an @var{n}-by-1 vector where @code{@var{n} = rows (@var{A})}\n\
00246 \n\
00247 @item disp\n\
00248 The level of diagnostic printout (0|1|2).  If @code{disp} is 0 then\n\
00249 diagnostics are disabled.  The default value is 0.\n\
00250 \n\
00251 @item cholB\n\
00252 Flag if @code{chol (@var{B})} is passed rather than @var{B}.  The default is\n\
00253 false.\n\
00254 \n\
00255 @item permB\n\
00256 The permutation vector of the Cholesky@tie{}factorization of @var{B} if\n\
00257 @code{cholB} is true.  That is @code{chol (@var{B}(permB, permB))}.  The\n\
00258 default is @code{1:@var{n}}.\n\
00259 \n\
00260 @end table\n\
00261 It is also possible to represent @var{A} by a function denoted @var{af}.\n\
00262 @var{af} must be followed by a scalar argument @var{n} defining the length\n\
00263 of the vector argument accepted by @var{af}.  @var{af} can be\n\
00264 a function handle, an inline function, or a string.  When @var{af} is a\n\
00265 string it holds the name of the function to use.\n\
00266 \n\
00267 @var{af} is a function of the form @code{y = af (x)}\n\
00268 where the required return value of @var{af} is determined by\n\
00269 the value of @var{sigma}.  The four possible forms are\n\
00270 \n\
00271 @table @code\n\
00272 @item A * x\n\
00273 if @var{sigma} is not given or is a string other than 'sm'.\n\
00274 \n\
00275 @item A \\ x\n\
00276 if @var{sigma} is 0 or 'sm'.\n\
00277 \n\
00278 @item (A - sigma * I) \\ x\n\
00279 for the standard eigenvalue problem, where @code{I} is the identity matrix of\n\
00280 the same size as @var{A}.\n\
00281 \n\
00282 @item (A - sigma * B) \\ x\n\
00283 for the general eigenvalue problem.\n\
00284 @end table\n\
00285 \n\
00286 The return arguments of @code{eigs} depend on the number of return arguments\n\
00287 requested.  With a single return argument, a vector @var{d} of length @var{k}\n\
00288 is returned containing the @var{k} eigenvalues that have been found.  With\n\
00289 two return arguments, @var{V} is a @var{n}-by-@var{k} matrix whose columns\n\
00290 are the @var{k} eigenvectors corresponding to the returned eigenvalues.  The\n\
00291 eigenvalues themselves are returned in @var{d} in the form of a\n\
00292 @var{n}-by-@var{k} matrix, where the elements on the diagonal are the\n\
00293 eigenvalues.\n\
00294 \n\
00295 Given a third return argument @var{flag}, @code{eigs} returns the status\n\
00296 of the convergence.  If @var{flag} is 0 then all eigenvalues have converged.\n\
00297 Any other value indicates a failure to converge.\n\
00298 \n\
00299 This function is based on the @sc{arpack} package, written by R. Lehoucq,\n\
00300 K. Maschhoff, D. Sorensen, and C. Yang.  For more information see\n\
00301 @url{http://www.caam.rice.edu/software/ARPACK/}.\n\
00302 \n\
00303 @seealso{eig, svds}\n\
00304 @end deftypefn")
00305 {
00306   octave_value_list retval;
00307 #ifdef HAVE_ARPACK
00308   int nargin = args.length ();
00309   std::string fcn_name;
00310   octave_idx_type n = 0;
00311   octave_idx_type k = 6;
00312   Complex sigma = 0.;
00313   double sigmar, sigmai;
00314   bool have_sigma = false;
00315   std::string typ = "LM";
00316   Matrix amm, bmm, bmt;
00317   ComplexMatrix acm, bcm, bct;
00318   SparseMatrix asmm, bsmm, bsmt;
00319   SparseComplexMatrix ascm, bscm, bsct;
00320   int b_arg = 0;
00321   bool have_b = false;
00322   bool have_a_fun = false;
00323   bool a_is_complex = false;
00324   bool b_is_complex = false;
00325   bool symmetric = false;
00326   bool sym_tested = false;
00327   bool cholB = false;
00328   bool a_is_sparse = false;
00329   ColumnVector permB;
00330   int arg_offset = 0;
00331   double tol = DBL_EPSILON;
00332   int maxit = 300;
00333   int disp = 0;
00334   octave_idx_type p = -1;
00335   ColumnVector resid;
00336   ComplexColumnVector cresid;
00337   octave_idx_type info = 1;
00338 
00339   warned_imaginary = false;
00340 
00341   unwind_protect frame;
00342 
00343   frame.protect_var (call_depth);
00344   call_depth++;
00345 
00346   if (call_depth > 1)
00347     {
00348       error ("eigs: invalid recursive call");
00349       if (fcn_name.length())
00350         clear_function (fcn_name);
00351       return retval;
00352     }
00353 
00354   if (nargin == 0)
00355     print_usage ();
00356   else if (args(0).is_function_handle () || args(0).is_inline_function ()
00357            || args(0).is_string ())
00358     {
00359       if (args(0).is_string ())
00360         {
00361           std::string name = args(0).string_value ();
00362           std::string fname = "function y = ";
00363           fcn_name = unique_symbol_name ("__eigs_fcn_");
00364           fname.append (fcn_name);
00365           fname.append ("(x) y = ");
00366           eigs_fcn = extract_function (args(0), "eigs", fcn_name, fname,
00367                                        "; endfunction");
00368         }
00369       else
00370         eigs_fcn = args(0).function_value ();
00371 
00372       if (!eigs_fcn)
00373         {
00374           error ("eigs: unknown function");
00375           return retval;
00376         }
00377 
00378       if (nargin < 2)
00379         {
00380           error ("eigs: incorrect number of arguments");
00381           return retval;
00382         }
00383       else
00384         {
00385           n = args(1).nint_value ();
00386           arg_offset = 1;
00387           have_a_fun = true;
00388         }
00389     }
00390   else
00391     {
00392       if (args(0).is_complex_type ())
00393         {
00394           if (args(0).is_sparse_type ())
00395             {
00396               ascm = (args(0).sparse_complex_matrix_value());
00397               a_is_sparse = true;
00398             }
00399           else
00400             acm = (args(0).complex_matrix_value());
00401           a_is_complex = true;
00402           symmetric = false; // ARPACK doesn't special case complex symmetric
00403           sym_tested = true;
00404         }
00405       else
00406         {
00407           if (args(0).is_sparse_type ())
00408             {
00409               asmm = (args(0).sparse_matrix_value());
00410               a_is_sparse = true;
00411             }
00412           else
00413             {
00414               amm = (args(0).matrix_value());
00415             }
00416         }
00417 
00418     }
00419 
00420   // Note hold off reading B till later to avoid issues of double
00421   // copies of the matrix if B is full/real while A is complex.
00422   if (!error_state && nargin > 1 + arg_offset &&
00423       !(args(1 + arg_offset).is_real_scalar ()))
00424     {
00425       if (args(1+arg_offset).is_complex_type ())
00426         {
00427           b_arg = 1+arg_offset;
00428           have_b = true;
00429           b_is_complex = true;
00430           arg_offset++;
00431         }
00432       else
00433         {
00434           b_arg = 1+arg_offset;
00435           have_b = true;
00436           arg_offset++;
00437         }
00438     }
00439 
00440   if (!error_state && nargin > (1+arg_offset))
00441     k = args(1+arg_offset).nint_value ();
00442 
00443   if (!error_state && nargin > (2+arg_offset))
00444     {
00445       if (args(2+arg_offset).is_string ())
00446         {
00447           typ = args(2+arg_offset).string_value ();
00448 
00449           // Use STL function to convert to upper case
00450           transform (typ.begin (), typ.end (), typ.begin (), toupper);
00451 
00452           sigma = 0.;
00453         }
00454       else
00455         {
00456           sigma = args(2+arg_offset).complex_value ();
00457 
00458           if (! error_state)
00459             have_sigma = true;
00460           else
00461             {
00462               error ("eigs: SIGMA must be a scalar or a string");
00463               return retval;
00464             }
00465         }
00466     }
00467 
00468   sigmar = std::real (sigma);
00469   sigmai = std::imag (sigma);
00470 
00471   if (!error_state && nargin > (3+arg_offset))
00472     {
00473       if (args(3+arg_offset).is_map ())
00474         {
00475           octave_scalar_map map = args(3+arg_offset).scalar_map_value ();
00476 
00477           if (! error_state)
00478             {
00479               octave_value tmp;
00480 
00481               // issym is ignored for complex matrix inputs
00482               tmp = map.getfield ("issym");
00483               if (tmp.is_defined () && !sym_tested)
00484                 {
00485                   symmetric = tmp.double_value () != 0.;
00486                   sym_tested = true;
00487                 }
00488 
00489               // isreal is ignored if A is not a function
00490               tmp = map.getfield ("isreal");
00491               if (tmp.is_defined () && have_a_fun)
00492                 a_is_complex = ! (tmp.double_value () != 0.);
00493 
00494               tmp = map.getfield ("tol");
00495               if (tmp.is_defined ())
00496                 tol = tmp.double_value ();
00497 
00498               tmp = map.getfield ("maxit");
00499               if (tmp.is_defined ())
00500                 maxit = tmp.nint_value ();
00501 
00502               tmp = map.getfield ("p");
00503               if (tmp.is_defined ())
00504                 p = tmp.nint_value ();
00505 
00506               tmp = map.getfield ("v0");
00507               if (tmp.is_defined ())
00508                 {
00509                   if (a_is_complex || b_is_complex)
00510                     cresid = ComplexColumnVector (tmp.complex_vector_value ());
00511                   else
00512                     resid = ColumnVector (tmp.vector_value ());
00513                 }
00514 
00515               tmp = map.getfield ("disp");
00516               if (tmp.is_defined ())
00517                 disp = tmp.nint_value ();
00518 
00519               tmp = map.getfield ("cholB");
00520               if (tmp.is_defined ())
00521                 cholB = tmp.double_value () != 0.;
00522 
00523               tmp = map.getfield ("permB");
00524               if (tmp.is_defined ())
00525                 permB = ColumnVector (tmp.vector_value ()) - 1.0;
00526             }
00527           else
00528             {
00529               error ("eigs: OPTS argument must be a scalar structure");
00530               return retval;
00531             }
00532         }
00533       else
00534         {
00535           error ("eigs: OPTS argument must be a structure");
00536           return retval;
00537         }
00538     }
00539 
00540   if (nargin > (4+arg_offset))
00541     {
00542       error ("eigs: incorrect number of arguments");
00543       return retval;
00544     }
00545 
00546   // Test undeclared (no issym) matrix inputs for symmetry
00547   if (!sym_tested && !have_a_fun)
00548     {
00549       if (a_is_sparse)
00550         symmetric = asmm.is_symmetric();
00551       else
00552         symmetric = amm.is_symmetric();
00553     }
00554 
00555   if (have_b)
00556     {
00557       if (a_is_complex || b_is_complex)
00558         {
00559           if (a_is_sparse)
00560             bscm = args(b_arg).sparse_complex_matrix_value ();
00561           else
00562             bcm = args(b_arg).complex_matrix_value ();
00563         }
00564       else
00565         {
00566           if (a_is_sparse)
00567             bsmm = args(b_arg).sparse_matrix_value ();
00568           else
00569             bmm = args(b_arg).matrix_value ();
00570         }
00571     }
00572 
00573   // Mode 1 for SM mode seems unstable for some reason.
00574   // Use Mode 3 instead, with sigma = 0.
00575   if (!error_state && !have_sigma && typ == "SM")
00576     have_sigma = true;
00577 
00578   if (!error_state)
00579     {
00580       octave_idx_type nconv;
00581       if (a_is_complex || b_is_complex)
00582         {
00583           ComplexMatrix eig_vec;
00584           ComplexColumnVector eig_val;
00585 
00586 
00587           if (have_a_fun)
00588             nconv = EigsComplexNonSymmetricFunc
00589               (eigs_complex_func, n, typ, sigma, k, p, info, eig_vec, eig_val,
00590                cresid, octave_stdout, tol, (nargout > 1), cholB, disp, maxit);
00591           else if (have_sigma)
00592             {
00593               if (a_is_sparse)
00594                 nconv = EigsComplexNonSymmetricMatrixShift
00595                   (ascm, sigma, k, p, info, eig_vec, eig_val, bscm, permB,
00596                    cresid, octave_stdout, tol, (nargout > 1), cholB, disp,
00597                    maxit);
00598               else
00599                 nconv = EigsComplexNonSymmetricMatrixShift
00600                   (acm, sigma, k, p, info, eig_vec, eig_val, bcm, permB, cresid,
00601                    octave_stdout, tol, (nargout > 1), cholB, disp, maxit);
00602             }
00603           else
00604             {
00605               if (a_is_sparse)
00606                 nconv = EigsComplexNonSymmetricMatrix
00607                   (ascm, typ, k, p, info, eig_vec, eig_val, bscm, permB, cresid,
00608                    octave_stdout, tol, (nargout > 1), cholB, disp, maxit);
00609               else
00610                 nconv = EigsComplexNonSymmetricMatrix
00611                   (acm, typ, k, p, info, eig_vec, eig_val, bcm, permB, cresid,
00612                    octave_stdout, tol, (nargout > 1), cholB, disp, maxit);
00613             }
00614 
00615           if (nargout < 2)
00616             retval (0) = eig_val;
00617           else
00618             {
00619               retval(2) = double (info);
00620               retval(1) = ComplexDiagMatrix (eig_val);
00621               retval(0) = eig_vec;
00622             }
00623         }
00624       else if (sigmai != 0.)
00625         {
00626           // Promote real problem to a complex one.
00627           ComplexMatrix eig_vec;
00628           ComplexColumnVector eig_val;
00629 
00630           if (have_a_fun)
00631             nconv = EigsComplexNonSymmetricFunc
00632               (eigs_complex_func, n, typ,  sigma, k, p, info, eig_vec, eig_val,
00633                cresid, octave_stdout, tol, (nargout > 1), cholB, disp, maxit);
00634           else
00635             {
00636               if (a_is_sparse)
00637                 nconv = EigsComplexNonSymmetricMatrixShift
00638                   (SparseComplexMatrix (asmm), sigma, k, p, info, eig_vec,
00639                    eig_val, SparseComplexMatrix (bsmm), permB, cresid,
00640                    octave_stdout, tol, (nargout > 1), cholB, disp, maxit);
00641               else
00642                 nconv = EigsComplexNonSymmetricMatrixShift
00643                   (ComplexMatrix (amm), sigma, k, p, info, eig_vec,
00644                    eig_val, ComplexMatrix (bmm), permB, cresid,
00645                    octave_stdout, tol, (nargout > 1), cholB, disp, maxit);
00646             }
00647 
00648           if (nargout < 2)
00649             retval (0) = eig_val;
00650           else
00651             {
00652               retval(2) = double (info);
00653               retval(1) = ComplexDiagMatrix (eig_val);
00654               retval(0) = eig_vec;
00655             }
00656         }
00657       else
00658         {
00659           if (symmetric)
00660             {
00661               Matrix eig_vec;
00662               ColumnVector eig_val;
00663 
00664               if (have_a_fun)
00665                 nconv = EigsRealSymmetricFunc
00666                   (eigs_func, n, typ, sigmar, k, p, info, eig_vec, eig_val,
00667                    resid, octave_stdout, tol, (nargout > 1), cholB, disp,
00668                    maxit);
00669               else if (have_sigma)
00670                 {
00671                   if (a_is_sparse)
00672                     nconv = EigsRealSymmetricMatrixShift
00673                       (asmm, sigmar, k, p, info, eig_vec, eig_val, bsmm, permB,
00674                        resid, octave_stdout, tol, (nargout > 1), cholB, disp,
00675                        maxit);
00676                   else
00677                     nconv = EigsRealSymmetricMatrixShift
00678                       (amm, sigmar, k, p, info, eig_vec, eig_val, bmm, permB,
00679                        resid, octave_stdout, tol, (nargout > 1), cholB, disp,
00680                        maxit);
00681                 }
00682               else
00683                 {
00684                   if (a_is_sparse)
00685                     nconv = EigsRealSymmetricMatrix
00686                       (asmm, typ, k, p, info, eig_vec, eig_val, bsmm, permB,
00687                        resid, octave_stdout, tol, (nargout > 1), cholB, disp,
00688                        maxit);
00689                   else
00690                     nconv = EigsRealSymmetricMatrix
00691                       (amm, typ, k, p, info, eig_vec, eig_val, bmm, permB,
00692                        resid, octave_stdout, tol, (nargout > 1), cholB, disp,
00693                        maxit);
00694                 }
00695 
00696               if (nargout < 2)
00697                 retval (0) = eig_val;
00698               else
00699                 {
00700                   retval(2) = double (info);
00701                   retval(1) = DiagMatrix (eig_val);
00702                   retval(0) = eig_vec;
00703                 }
00704             }
00705           else
00706             {
00707               ComplexMatrix eig_vec;
00708               ComplexColumnVector eig_val;
00709 
00710               if (have_a_fun)
00711                 nconv = EigsRealNonSymmetricFunc
00712                   (eigs_func, n, typ, sigmar, k, p, info, eig_vec, eig_val,
00713                    resid, octave_stdout, tol, (nargout > 1), cholB, disp,
00714                    maxit);
00715               else if (have_sigma)
00716                 {
00717                   if (a_is_sparse)
00718                     nconv = EigsRealNonSymmetricMatrixShift
00719                       (asmm, sigmar, k, p, info, eig_vec, eig_val, bsmm, permB,
00720                        resid, octave_stdout, tol, (nargout > 1), cholB, disp,
00721                        maxit);
00722                   else
00723                     nconv = EigsRealNonSymmetricMatrixShift
00724                       (amm, sigmar, k, p, info, eig_vec, eig_val, bmm, permB,
00725                        resid, octave_stdout, tol, (nargout > 1), cholB, disp,
00726                        maxit);
00727                 }
00728               else
00729                 {
00730                   if (a_is_sparse)
00731                     nconv = EigsRealNonSymmetricMatrix
00732                       (asmm, typ, k, p, info, eig_vec, eig_val, bsmm, permB,
00733                        resid, octave_stdout, tol, (nargout > 1), cholB, disp,
00734                        maxit);
00735                   else
00736                     nconv = EigsRealNonSymmetricMatrix
00737                       (amm, typ, k, p, info, eig_vec, eig_val, bmm, permB,
00738                        resid, octave_stdout, tol, (nargout > 1), cholB, disp,
00739                        maxit);
00740                 }
00741 
00742               if (nargout < 2)
00743                 retval (0) = eig_val;
00744               else
00745                 {
00746                   retval(2) = double (info);
00747                   retval(1) = ComplexDiagMatrix (eig_val);
00748                   retval(0) = eig_vec;
00749                 }
00750             }
00751         }
00752 
00753       if (nconv <= 0)
00754         warning ("eigs: None of the %d requested eigenvalues converged", k);
00755       else if (nconv < k)
00756         warning ("eigs: Only %d of the %d requested eigenvalues converged",
00757                  nconv, k);
00758     }
00759 
00760   if (! fcn_name.empty ())
00761     clear_function (fcn_name);
00762 #else
00763   error ("eigs: not available in this version of Octave");
00764 #endif
00765 
00766   return retval;
00767 }
00768 
00769 /* #### SPARSE MATRIX VERSIONS #### */
00770 
00771 /*
00772 
00773 %% Real positive definite tests, n must be even
00774 %!shared n, k, A, d0, d2
00775 %! n = 20;
00776 %! k = 4;
00777 %! A = sparse([3:n,1:n,1:(n-2)],[1:(n-2),1:n,3:n],[ones(1,n-2),4*ones(1,n),ones(1,n-2)]);
00778 %! d0 = eig (A);
00779 %! d2 = sort (d0);
00780 %! [~, idx] = sort (abs(d0));
00781 %! d0 = d0(idx);
00782 %! rand("state", 42); % initialize generator to make eigs behavior reproducible
00783 %!testif HAVE_ARPACK
00784 %! d1 = eigs (A, k);
00785 %! assert (d1, d0(end:-1:(end-k+1)), 1e-11);
00786 %!testif HAVE_ARPACK
00787 %! d1 = eigs (A,k+1);
00788 %! assert (d1, d0(end:-1:(end-k)),1e-11);
00789 %!testif HAVE_ARPACK
00790 %! d1 = eigs (A, k, 'lm');
00791 %! assert (d1, d0(end:-1:(end-k+1)), 1e-11);
00792 %!testif HAVE_ARPACK, HAVE_UMFPACK
00793 %! d1 = eigs (A, k, 'sm');
00794 %! assert (d1, d0(k:-1:1), 1e-11);
00795 %!testif HAVE_ARPACK
00796 %! d1 = eigs (A, k, 'la');
00797 %! assert (d1, d2(end:-1:(end-k+1)), 1e-11);
00798 %!testif HAVE_ARPACK
00799 %! d1 = eigs (A, k, 'sa');
00800 %! assert (d1, d2(1:k), 1e-11);
00801 %!testif HAVE_ARPACK
00802 %! d1 = eigs (A, k, 'be');
00803 %! assert (d1, d2([1:floor(k/2), (end - ceil(k/2) + 1):end]), 1e-11);
00804 %!testif HAVE_ARPACK
00805 %! d1 = eigs (A, k+1, 'be');
00806 %! assert (d1, d2([1:floor((k+1)/2), (end - ceil((k+1)/2) + 1):end]), 1e-11);
00807 %!testif HAVE_ARPACK, HAVE_UMFPACK
00808 %! d1 = eigs (A, k, 4.1);
00809 %! [~,idx0] = sort (abs(d0 - 4.1));
00810 %! [~,idx1] = sort (abs(d1 - 4.1));
00811 %! assert (d1(idx1), d0(idx0(1:k)), 1e-11);
00812 %!testif HAVE_ARPACK, HAVE_CHOLMOD
00813 %! d1 = eigs(A, speye(n), k, 'lm');
00814 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
00815 %!testif HAVE_ARPACK, HAVE_UMFPACK
00816 %! assert (eigs(A,k,4.1), eigs(A,speye(n),k,4.1), 1e-11);
00817 %!testif HAVE_ARPACK
00818 %! opts.cholB=true;
00819 %! d1 = eigs(A, speye(n), k, 'lm', opts);
00820 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
00821 %!testif HAVE_ARPACK
00822 %! opts.cholB=true;
00823 %! q = [2:n,1];
00824 %! opts.permB=q;
00825 %! d1 = eigs(A, speye(n)(q,q), k, 'lm', opts);
00826 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
00827 %!testif HAVE_ARPACK, HAVE_UMFPACK
00828 %! opts.cholB=true;
00829 %! d1 = eigs(A, speye(n), k, 4.1, opts);
00830 %! assert (abs(d1), eigs(A,k,4.1), 1e-11);
00831 %!testif HAVE_ARPACK, HAVE_UMFPACK
00832 %! opts.cholB=true;
00833 %! q = [2:n,1];
00834 %! opts.permB=q;
00835 %! d1 = eigs(A, speye(n)(q,q), k, 4.1, opts);
00836 %! assert (abs(d1), eigs(A,k,4.1), 1e-11);
00837 %!testif HAVE_ARPACK, HAVE_UMFPACK
00838 %! assert (eigs(A,k,4.1), eigs(A,speye(n),k,4.1), 1e-11);
00839 %!testif HAVE_ARPACK
00840 %! fn = @(x) A * x;
00841 %! opts.issym = 1; opts.isreal = 1;
00842 %! d1 = eigs (fn, n, k, 'lm', opts);
00843 %! assert (d1, d0(end:-1:(end-k+1)), 1e-11);
00844 %!testif HAVE_ARPACK
00845 %! fn = @(x) A \ x;
00846 %! opts.issym = 1; opts.isreal = 1;
00847 %! d1 = eigs (fn, n, k, 'sm', opts);
00848 %! assert (d1, d0(k:-1:1), 1e-11);
00849 %!testif HAVE_ARPACK, HAVE_UMFPACK
00850 %! fn = @(x) (A - 4.1 * eye(n)) \ x;
00851 %! opts.issym = 1; opts.isreal = 1;
00852 %! d1 = eigs (fn, n, k, 4.1, opts);
00853 %! assert (d1, eigs(A,k,4.1), 1e-11);
00854 %!testif HAVE_ARPACK
00855 %! AA = speye (10);
00856 %! fn = @(x) AA * x;
00857 %! opts.issym = 1; opts.isreal = 1;
00858 %! assert (eigs (fn, 10, AA, 3, 'lm', opts), [1; 1; 1],10*eps);
00859 %!testif HAVE_ARPACK
00860 %! [v1,d1] = eigs(A, k, 'lm');
00861 %! d1 = diag(d1);
00862 %! for i=1:k
00863 %!  assert(max(abs((A - d1(i)*speye(n))*v1(:,i))),0.,1e-11)
00864 %! endfor
00865 %!testif HAVE_ARPACK, HAVE_UMFPACK
00866 %! [v1,d1] = eigs(A, k, 'sm');
00867 %! d1 = diag(d1);
00868 %! for i=1:k
00869 %!  assert(max(abs((A - d1(i)*speye(n))*v1(:,i))),0.,1e-11)
00870 %! endfor
00871 %!testif HAVE_ARPACK
00872 %! [v1,d1] = eigs(A, k, 'la');
00873 %! d1 = diag(d1);
00874 %! for i=1:k
00875 %!  assert(max(abs((A - d1(i)*speye(n))*v1(:,i))),0.,1e-11)
00876 %! endfor
00877 %!testif HAVE_ARPACK
00878 %! [v1,d1] = eigs(A, k, 'sa');
00879 %! d1 = diag(d1);
00880 %! for i=1:k
00881 %!  assert(max(abs((A - d1(i)*speye(n))*v1(:,i))),0.,1e-11)
00882 %! endfor
00883 %!testif HAVE_ARPACK
00884 %! [v1,d1] = eigs(A, k, 'be');
00885 %! d1 = diag(d1);
00886 %! for i=1:k
00887 %!  assert(max(abs((A - d1(i)*speye(n))*v1(:,i))),0.,1e-11)
00888 %! endfor
00889 
00890 */
00891 
00892 /*
00893 
00894 %% Real unsymmetric tests
00895 %!shared n, k, A, d0
00896 %! n = 20;
00897 %! k = 4;
00898 %! A =  sparse([3:n,1:n,1:(n-2)],[1:(n-2),1:n,3:n],[ones(1,n-2),1:n,-ones(1,n-2)]);
00899 %! d0 = eig (A);
00900 %! [~, idx] = sort (abs(d0));
00901 %! d0 = d0(idx);
00902 %! rand("state", 42); % initialize generator to make eigs behavior reproducible
00903 %!testif HAVE_ARPACK
00904 %! d1 = eigs (A, k);
00905 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
00906 %!testif HAVE_ARPACK
00907 %! d1 = eigs (A,k+1);
00908 %! assert (abs(d1), abs(d0(end:-1:(end-k))),1e-11);
00909 %!testif HAVE_ARPACK
00910 %! d1 = eigs (A, k, 'lm');
00911 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
00912 %!testif HAVE_ARPACK, HAVE_UMFPACK
00913 %! d1 = eigs (A, k, 'sm');
00914 %! assert (abs(d1), abs(d0(1:k)), 1e-11);
00915 %!testif HAVE_ARPACK
00916 %! d1 = eigs (A, k, 'lr');
00917 %! [~, idx] = sort (real(d0));
00918 %! d2 = d0(idx);
00919 %! assert (real(d1), real(d2(end:-1:(end-k+1))), 1e-11);
00920 %!testif HAVE_ARPACK
00921 %! d1 = eigs (A, k, 'sr');
00922 %! [~, idx] = sort (real(abs(d0)));
00923 %! d2 = d0(idx);
00924 %! assert (real(d1), real(d2(1:k)), 1e-11);
00925 %!testif HAVE_ARPACK
00926 %! d1 = eigs (A, k, 'li');
00927 %! [~, idx] = sort (imag(abs(d0)));
00928 %! d2 = d0(idx);
00929 %! assert (sort(imag(d1)), sort(imag(d2(end:-1:(end-k+1)))), 1e-11);
00930 %!testif HAVE_ARPACK
00931 %! d1 = eigs (A, k, 'si');
00932 %! [~, idx] = sort (imag(abs(d0)));
00933 %! d2 = d0(idx);
00934 %! assert (sort(imag(d1)), sort(imag(d2(1:k))), 1e-11);
00935 %!testif HAVE_ARPACK, HAVE_UMFPACK
00936 %! d1 = eigs (A, k, 4.1);
00937 %! [~,idx0] = sort (abs(d0 - 4.1));
00938 %! [~,idx1] = sort (abs(d1 - 4.1));
00939 %! assert (abs(d1(idx1)), abs(d0(idx0(1:k))), 1e-11);
00940 %! assert (sort(imag(d1(idx1))), sort(imag(d0(idx0(1:k)))), 1e-11);
00941 %!testif HAVE_ARPACK, HAVE_CHOLMOD
00942 %! d1 = eigs(A, speye(n), k, 'lm');
00943 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
00944 %!testif HAVE_ARPACK
00945 %! opts.cholB=true;
00946 %! d1 = eigs(A, speye(n), k, 'lm', opts);
00947 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
00948 %!testif HAVE_ARPACK
00949 %! opts.cholB=true;
00950 %! q = [2:n,1];
00951 %! opts.permB=q;
00952 %! d1 = eigs(A, speye(n)(q,q), k, 'lm', opts);
00953 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
00954 %!testif HAVE_ARPACK, HAVE_UMFPACK
00955 %! opts.cholB=true;
00956 %! d1 = eigs(A, speye(n), k, 4.1, opts);
00957 %! assert (abs(d1), eigs(A,k,4.1), 1e-11);
00958 %!testif HAVE_ARPACK, HAVE_UMFPACK
00959 %! opts.cholB=true;
00960 %! q = [2:n,1];
00961 %! opts.permB=q;
00962 %! d1 = eigs(A, speye(n)(q,q), k, 4.1, opts);
00963 %! assert (abs(d1), eigs(A,k,4.1), 1e-11);
00964 %!testif HAVE_ARPACK, HAVE_UMFPACK
00965 %! assert (abs(eigs(A,k,4.1)), abs(eigs(A,speye(n),k,4.1)), 1e-11);
00966 %!testif HAVE_ARPACK, HAVE_UMFPACK
00967 %! assert (sort(imag(eigs(A,k,4.1))), sort(imag(eigs(A,speye(n),k,4.1))), 1e-11);
00968 %!testif HAVE_ARPACK
00969 %! fn = @(x) A * x;
00970 %! opts.issym = 0; opts.isreal = 1;
00971 %! d1 = eigs (fn, n, k, 'lm', opts);
00972 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
00973 %!testif HAVE_ARPACK
00974 %! fn = @(x) A \ x;
00975 %! opts.issym = 0; opts.isreal = 1;
00976 %! d1 = eigs (fn, n, k, 'sm', opts);
00977 %! assert (abs(d1), d0(1:k), 1e-11);
00978 %!testif HAVE_ARPACK, HAVE_UMFPACK
00979 %! fn = @(x) (A - 4.1 * eye(n)) \ x;
00980 %! opts.issym = 0; opts.isreal = 1;
00981 %! d1 = eigs (fn, n, k, 4.1, opts);
00982 %! assert (abs(d1), eigs(A,k,4.1), 1e-11);
00983 %!testif HAVE_ARPACK
00984 %! [v1,d1] = eigs(A, k, 'lm');
00985 %! d1 = diag(d1);
00986 %! for i=1:k
00987 %!  assert(max(abs((A - d1(i)*speye(n))*v1(:,i))),0.,1e-11)
00988 %! endfor
00989 %!testif HAVE_ARPACK, HAVE_UMFPACK
00990 %! [v1,d1] = eigs(A, k, 'sm');
00991 %! d1 = diag(d1);
00992 %! for i=1:k
00993 %!  assert(max(abs((A - d1(i)*speye(n))*v1(:,i))),0.,1e-11)
00994 %! endfor
00995 %!testif HAVE_ARPACK
00996 %! [v1,d1] = eigs(A, k, 'lr');
00997 %! d1 = diag(d1);
00998 %! for i=1:k
00999 %!  assert(max(abs((A - d1(i)*speye(n))*v1(:,i))),0.,1e-11)
01000 %! endfor
01001 %!testif HAVE_ARPACK
01002 %! [v1,d1] = eigs(A, k, 'sr');
01003 %! d1 = diag(d1);
01004 %! for i=1:k
01005 %!  assert(max(abs((A - d1(i)*speye(n))*v1(:,i))),0.,1e-11)
01006 %! endfor
01007 %!testif HAVE_ARPACK
01008 %! [v1,d1] = eigs(A, k, 'li');
01009 %! d1 = diag(d1);
01010 %! for i=1:k
01011 %!  assert(max(abs((A - d1(i)*speye(n))*v1(:,i))),0.,1e-11)
01012 %! endfor
01013 %!testif HAVE_ARPACK
01014 %! [v1,d1] = eigs(A, k, 'si');
01015 %! d1 = diag(d1);
01016 %! for i=1:k
01017 %!  assert(max(abs((A - d1(i)*speye(n))*v1(:,i))),0.,1e-11)
01018 %! endfor
01019 
01020 */
01021 
01022 /*
01023 
01024 %% Complex hermitian tests
01025 %!shared n, k, A, d0
01026 %! n = 20;
01027 %! k = 4;
01028 %! A = sparse([3:n,1:n,1:(n-2)],[1:(n-2),1:n,3:n],[1i*ones(1,n-2),4*ones(1,n),-1i*ones(1,n-2)]);
01029 %! d0 = eig (A);
01030 %! [~, idx] = sort (abs(d0));
01031 %! d0 = d0(idx);
01032 %! rand("state", 42); % initialize generator to make eigs behavior reproducible
01033 %!testif HAVE_ARPACK
01034 %! d1 = eigs (A, k);
01035 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01036 %!testif HAVE_ARPACK
01037 %! d1 = eigs (A,k+1);
01038 %! assert (abs(d1), abs(d0(end:-1:(end-k))),1e-11);
01039 %!testif HAVE_ARPACK
01040 %! d1 = eigs (A, k, 'lm');
01041 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01042 %!testif HAVE_ARPACK, HAVE_UMFPACK
01043 %! d1 = eigs (A, k, 'sm');
01044 %! assert (abs(d1), abs(d0(1:k)), 1e-11);
01045 %!testif HAVE_ARPACK
01046 %! d1 = eigs (A, k, 'lr');
01047 %! [~, idx] = sort (real(abs(d0)));
01048 %! d2 = d0(idx);
01049 %! assert (real(d1), real(d2(end:-1:(end-k+1))), 1e-11);
01050 %!testif HAVE_ARPACK
01051 %! d1 = eigs (A, k, 'sr');
01052 %! [~, idx] = sort (real(abs(d0)));
01053 %! d2 = d0(idx);
01054 %! assert (real(d1), real(d2(1:k)), 1e-11);
01055 %!testif HAVE_ARPACK
01056 %! d1 = eigs (A, k, 'li');
01057 %! [~, idx] = sort (imag(abs(d0)));
01058 %! d2 = d0(idx);
01059 %! assert (sort(imag(d1)), sort(imag(d2(end:-1:(end-k+1)))), 1e-11);
01060 %!testif HAVE_ARPACK
01061 %! d1 = eigs (A, k, 'si');
01062 %! [~, idx] = sort (imag(abs(d0)));
01063 %! d2 = d0(idx);
01064 %! assert (sort(imag(d1)), sort(imag(d2(1:k))), 1e-11);
01065 %!testif HAVE_ARPACK, HAVE_UMFPACK
01066 %! d1 = eigs (A, k, 4.1);
01067 %! [~,idx0] = sort (abs(d0 - 4.1));
01068 %! [~,idx1] = sort (abs(d1 - 4.1));
01069 %! assert (abs(d1(idx1)), abs(d0(idx0(1:k))), 1e-11);
01070 %! assert (sort(imag(d1(idx1))), sort(imag(d0(idx0(1:k)))), 1e-11);
01071 %!testif HAVE_ARPACK, HAVE_CHOLMOD
01072 %! d1 = eigs(A, speye(n), k, 'lm');
01073 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01074 %!testif HAVE_ARPACK
01075 %! opts.cholB=true;
01076 %! d1 = eigs(A, speye(n), k, 'lm', opts);
01077 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01078 %!testif HAVE_ARPACK
01079 %! opts.cholB=true;
01080 %! q = [2:n,1];
01081 %! opts.permB=q;
01082 %! d1 = eigs(A, speye(n)(q,q), k, 'lm', opts);
01083 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01084 %!testif HAVE_ARPACK, HAVE_UMFPACK
01085 %! opts.cholB=true;
01086 %! d1 = eigs(A, speye(n), k, 4.1, opts);
01087 %! assert (abs(abs(d1)), abs(eigs(A,k,4.1)), 1e-11);
01088 %! assert (sort(imag(abs(d1))), sort(imag(eigs(A,k,4.1))), 1e-11);
01089 %!testif HAVE_ARPACK, HAVE_UMFPACK
01090 %! opts.cholB=true;
01091 %! q = [2:n,1];
01092 %! opts.permB=q;
01093 %! d1 = eigs(A, speye(n)(q,q), k, 4.1, opts);
01094 %! assert (abs(abs(d1)), abs(eigs(A,k,4.1)), 1e-11);
01095 %! assert (sort(imag(abs(d1))), sort(imag(eigs(A,k,4.1))), 1e-11);
01096 %!testif HAVE_ARPACK, HAVE_UMFPACK
01097 %! assert (abs(eigs(A,k,4.1)), abs(eigs(A,speye(n),k,4.1)), 1e-11);
01098 %!testif HAVE_ARPACK, HAVE_UMFPACK
01099 %! assert (sort(imag(eigs(A,k,4.1))), sort(imag(eigs(A,speye(n),k,4.1))), 1e-11);
01100 %!testif HAVE_ARPACK
01101 %! fn = @(x) A * x;
01102 %! opts.issym = 0; opts.isreal = 0;
01103 %! d1 = eigs (fn, n, k, 'lm', opts);
01104 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01105 %!testif HAVE_ARPACK
01106 %! fn = @(x) A \ x;
01107 %! opts.issym = 0; opts.isreal = 0;
01108 %! d1 = eigs (fn, n, k, 'sm', opts);
01109 %! assert (abs(d1), d0(1:k), 1e-11);
01110 %!testif HAVE_ARPACK, HAVE_UMFPACK
01111 %! fn = @(x) (A - 4.1 * eye(n)) \ x;
01112 %! opts.issym = 0; opts.isreal = 0;
01113 %! d1 = eigs (fn, n, k, 4.1, opts);
01114 %! assert (abs(d1), eigs(A,k,4.1), 1e-11);
01115 %!testif HAVE_ARPACK
01116 %! [v1,d1] = eigs(A, k, 'lm');
01117 %! d1 = diag(d1);
01118 %! for i=1:k
01119 %!  assert(max(abs((A - d1(i)*speye(n))*v1(:,i))),0.,1e-11)
01120 %! endfor
01121 %!testif HAVE_ARPACK, HAVE_UMFPACK
01122 %! [v1,d1] = eigs(A, k, 'sm');
01123 %! d1 = diag(d1);
01124 %! for i=1:k
01125 %!  assert(max(abs((A - d1(i)*speye(n))*v1(:,i))),0.,1e-11)
01126 %! endfor
01127 %!testif HAVE_ARPACK
01128 %! [v1,d1] = eigs(A, k, 'lr');
01129 %! d1 = diag(d1);
01130 %! for i=1:k
01131 %!  assert(max(abs((A - d1(i)*speye(n))*v1(:,i))),0.,1e-11)
01132 %! endfor
01133 %!testif HAVE_ARPACK
01134 %! [v1,d1] = eigs(A, k, 'sr');
01135 %! d1 = diag(d1);
01136 %! for i=1:k
01137 %!  assert(max(abs((A - d1(i)*speye(n))*v1(:,i))),0.,1e-11)
01138 %! endfor
01139 %!testif HAVE_ARPACK
01140 %! [v1,d1] = eigs(A, k, 'li');
01141 %! d1 = diag(d1);
01142 %! for i=1:k
01143 %!  assert(max(abs((A - d1(i)*speye(n))*v1(:,i))),0.,1e-11)
01144 %! endfor
01145 %!testif HAVE_ARPACK
01146 %! [v1,d1] = eigs(A, k, 'si');
01147 %! d1 = diag(d1);
01148 %! for i=1:k
01149 %!  assert(max(abs((A - d1(i)*speye(n))*v1(:,i))),0.,1e-11)
01150 %! endfor
01151 
01152 */
01153 
01154 /* #### FULL MATRIX VERSIONS #### */
01155 
01156 /*
01157 
01158 %% Real positive definite tests, n must be even
01159 %!shared n, k, A, d0, d2
01160 %! n = 20;
01161 %! k = 4;
01162 %! A = full(sparse([3:n,1:n,1:(n-2)],[1:(n-2),1:n,3:n],[ones(1,n-2),4*ones(1,n),ones(1,n-2)]));
01163 %! d0 = eig (A);
01164 %! d2 = sort (d0);
01165 %! [~, idx] = sort (abs(d0));
01166 %! d0 = d0(idx);
01167 %! rand("state", 42); % initialize generator to make eigs behavior reproducible
01168 %!testif HAVE_ARPACK
01169 %! d1 = eigs (A, k);
01170 %! assert (d1, d0(end:-1:(end-k+1)), 1e-11);
01171 %!testif HAVE_ARPACK
01172 %! d1 = eigs (A,k+1);
01173 %! assert (d1, d0(end:-1:(end-k)),1e-11);
01174 %!testif HAVE_ARPACK
01175 %! d1 = eigs (A, k, 'lm');
01176 %! assert (d1, d0(end:-1:(end-k+1)), 1e-11);
01177 %!testif HAVE_ARPACK
01178 %! d1 = eigs (A, k, 'sm');
01179 %! assert (d1, d0(k:-1:1), 1e-11);
01180 %!testif HAVE_ARPACK
01181 %! d1 = eigs (A, k, 'la');
01182 %! assert (d1, d2(end:-1:(end-k+1)), 1e-11);
01183 %!testif HAVE_ARPACK
01184 %! d1 = eigs (A, k, 'sa');
01185 %! assert (d1, d2(1:k), 1e-11);
01186 %!testif HAVE_ARPACK
01187 %! d1 = eigs (A, k, 'be');
01188 %! assert (d1, d2([1:floor(k/2), (end - ceil(k/2) + 1):end]), 1e-11);
01189 %!testif HAVE_ARPACK
01190 %! d1 = eigs (A, k+1, 'be');
01191 %! assert (d1, d2([1:floor((k+1)/2), (end - ceil((k+1)/2) + 1):end]), 1e-11);
01192 %!testif HAVE_ARPACK
01193 %! d1 = eigs (A, k, 4.1);
01194 %! [~,idx0] = sort (abs(d0 - 4.1));
01195 %! [~,idx1] = sort (abs(d1 - 4.1));
01196 %! assert (d1(idx1), d0(idx0(1:k)), 1e-11);
01197 %!testif HAVE_ARPACK
01198 %! d1 = eigs(A, eye(n), k, 'lm');
01199 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01200 %!testif HAVE_ARPACK
01201 %! assert (eigs(A,k,4.1), eigs(A,eye(n),k,4.1), 1e-11);
01202 %!testif HAVE_ARPACK
01203 %! opts.cholB=true;
01204 %! d1 = eigs(A, eye(n), k, 'lm', opts);
01205 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01206 %!testif HAVE_ARPACK
01207 %! opts.cholB=true;
01208 %! q = [2:n,1];
01209 %! opts.permB=q;
01210 %! d1 = eigs(A, eye(n)(q,q), k, 'lm', opts);
01211 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01212 %!testif HAVE_ARPACK
01213 %! opts.cholB=true;
01214 %! d1 = eigs(A, eye(n), k, 4.1, opts);
01215 %! assert (abs(d1), eigs(A,k,4.1), 1e-11);
01216 %!testif HAVE_ARPACK
01217 %! opts.cholB=true;
01218 %! q = [2:n,1];
01219 %! opts.permB=q;
01220 %! d1 = eigs(A, eye(n)(q,q), k, 4.1, opts);
01221 %! assert (abs(d1), eigs(A,k,4.1), 1e-11);
01222 %!testif HAVE_ARPACK
01223 %! assert (eigs(A,k,4.1), eigs(A,eye(n),k,4.1), 1e-11);
01224 %!testif HAVE_ARPACK
01225 %! fn = @(x) A * x;
01226 %! opts.issym = 1; opts.isreal = 1;
01227 %! d1 = eigs (fn, n, k, 'lm', opts);
01228 %! assert (d1, d0(end:-1:(end-k+1)), 1e-11);
01229 %!testif HAVE_ARPACK
01230 %! fn = @(x) A \ x;
01231 %! opts.issym = 1; opts.isreal = 1;
01232 %! d1 = eigs (fn, n, k, 'sm', opts);
01233 %! assert (d1, d0(k:-1:1), 1e-11);
01234 %!testif HAVE_ARPACK
01235 %! fn = @(x) (A - 4.1 * eye(n)) \ x;
01236 %! opts.issym = 1; opts.isreal = 1;
01237 %! d1 = eigs (fn, n, k, 4.1, opts);
01238 %! assert (d1, eigs(A,k,4.1), 1e-11);
01239 %!testif HAVE_ARPACK
01240 %! [v1,d1] = eigs(A, k, 'lm');
01241 %! d1 = diag(d1);
01242 %! for i=1:k
01243 %!  assert(max(abs((A - d1(i)*eye(n))*v1(:,i))),0.,1e-11)
01244 %! endfor
01245 %!testif HAVE_ARPACK
01246 %! [v1,d1] = eigs(A, k, 'sm');
01247 %! d1 = diag(d1);
01248 %! for i=1:k
01249 %!  assert(max(abs((A - d1(i)*eye(n))*v1(:,i))),0.,1e-11)
01250 %! endfor
01251 %!testif HAVE_ARPACK
01252 %! [v1,d1] = eigs(A, k, 'la');
01253 %! d1 = diag(d1);
01254 %! for i=1:k
01255 %!  assert(max(abs((A - d1(i)*eye(n))*v1(:,i))),0.,1e-11)
01256 %! endfor
01257 %!testif HAVE_ARPACK
01258 %! [v1,d1] = eigs(A, k, 'sa');
01259 %! d1 = diag(d1);
01260 %! for i=1:k
01261 %!  assert(max(abs((A - d1(i)*eye(n))*v1(:,i))),0.,1e-11)
01262 %! endfor
01263 %!testif HAVE_ARPACK
01264 %! [v1,d1] = eigs(A, k, 'be');
01265 %! d1 = diag(d1);
01266 %! for i=1:k
01267 %!  assert(max(abs((A - d1(i)*eye(n))*v1(:,i))),0.,1e-11)
01268 %! endfor
01269 
01270 */
01271 
01272 /*
01273 
01274 %% Real unsymmetric tests
01275 %!shared n, k, A, d0
01276 %! n = 20;
01277 %! k = 4;
01278 %! A =  full(sparse([3:n,1:n,1:(n-2)],[1:(n-2),1:n,3:n],[ones(1,n-2),1:n,-ones(1,n-2)]));
01279 %! d0 = eig (A);
01280 %! [~, idx] = sort (abs(d0));
01281 %! d0 = d0(idx);
01282 %! rand("state", 42); % initialize generator to make eigs behavior reproducible
01283 %!testif HAVE_ARPACK
01284 %! d1 = eigs (A, k);
01285 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01286 %!testif HAVE_ARPACK
01287 %! d1 = eigs (A,k+1);
01288 %! assert (abs(d1), abs(d0(end:-1:(end-k))),1e-11);
01289 %!testif HAVE_ARPACK
01290 %! d1 = eigs (A, k, 'lm');
01291 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01292 %!testif HAVE_ARPACK
01293 %! d1 = eigs (A, k, 'sm');
01294 %! assert (abs(d1), abs(d0(1:k)), 1e-11);
01295 %!testif HAVE_ARPACK
01296 %! d1 = eigs (A, k, 'lr');
01297 %! [~, idx] = sort (real(d0));
01298 %! d2 = d0(idx);
01299 %! assert (real(d1), real(d2(end:-1:(end-k+1))), 1e-11);
01300 %!testif HAVE_ARPACK
01301 %! d1 = eigs (A, k, 'sr');
01302 %! [~, idx] = sort (real(abs(d0)));
01303 %! d2 = d0(idx);
01304 %! assert (real(d1), real(d2(1:k)), 1e-11);
01305 %!testif HAVE_ARPACK
01306 %! d1 = eigs (A, k, 'li');
01307 %! [~, idx] = sort (imag(abs(d0)));
01308 %! d2 = d0(idx);
01309 %! assert (sort(imag(d1)), sort(imag(d2(end:-1:(end-k+1)))), 1e-11);
01310 %!testif HAVE_ARPACK
01311 %! d1 = eigs (A, k, 'si');
01312 %! [~, idx] = sort (imag(abs(d0)));
01313 %! d2 = d0(idx);
01314 %! assert (sort(imag(d1)), sort(imag(d2(1:k))), 1e-11);
01315 %!testif HAVE_ARPACK
01316 %! d1 = eigs (A, k, 4.1);
01317 %! [~,idx0] = sort (abs(d0 - 4.1));
01318 %! [~,idx1] = sort (abs(d1 - 4.1));
01319 %! assert (abs(d1(idx1)), abs(d0(idx0(1:k))), 1e-11);
01320 %! assert (sort(imag(d1(idx1))), sort(imag(d0(idx0(1:k)))), 1e-11);
01321 %!testif HAVE_ARPACK
01322 %! d1 = eigs(A, eye(n), k, 'lm');
01323 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01324 %!testif HAVE_ARPACK
01325 %! opts.cholB=true;
01326 %! d1 = eigs(A, eye(n), k, 'lm', opts);
01327 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01328 %!testif HAVE_ARPACK
01329 %! opts.cholB=true;
01330 %! q = [2:n,1];
01331 %! opts.permB=q;
01332 %! d1 = eigs(A, eye(n)(q,q), k, 'lm', opts);
01333 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01334 %!testif HAVE_ARPACK
01335 %! opts.cholB=true;
01336 %! d1 = eigs(A, eye(n), k, 4.1, opts);
01337 %! assert (abs(d1), eigs(A,k,4.1), 1e-11);
01338 %!testif HAVE_ARPACK
01339 %! opts.cholB=true;
01340 %! q = [2:n,1];
01341 %! opts.permB=q;
01342 %! d1 = eigs(A, eye(n)(q,q), k, 4.1, opts);
01343 %! assert (abs(d1), eigs(A,k,4.1), 1e-11);
01344 %!testif HAVE_ARPACK
01345 %! assert (abs(eigs(A,k,4.1)), abs(eigs(A,eye(n),k,4.1)), 1e-11);
01346 %!testif HAVE_ARPACK
01347 %! assert (sort(imag(eigs(A,k,4.1))), sort(imag(eigs(A,eye(n),k,4.1))), 1e-11);
01348 %!testif HAVE_ARPACK
01349 %! fn = @(x) A * x;
01350 %! opts.issym = 0; opts.isreal = 1;
01351 %! d1 = eigs (fn, n, k, 'lm', opts);
01352 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01353 %!testif HAVE_ARPACK
01354 %! fn = @(x) A \ x;
01355 %! opts.issym = 0; opts.isreal = 1;
01356 %! d1 = eigs (fn, n, k, 'sm', opts);
01357 %! assert (abs(d1), d0(1:k), 1e-11);
01358 %!testif HAVE_ARPACK
01359 %! fn = @(x) (A - 4.1 * eye(n)) \ x;
01360 %! opts.issym = 0; opts.isreal = 1;
01361 %! d1 = eigs (fn, n, k, 4.1, opts);
01362 %! assert (abs(d1), eigs(A,k,4.1), 1e-11);
01363 %!testif HAVE_ARPACK
01364 %! [v1,d1] = eigs(A, k, 'lm');
01365 %! d1 = diag(d1);
01366 %! for i=1:k
01367 %!  assert(max(abs((A - d1(i)*eye(n))*v1(:,i))),0.,1e-11)
01368 %! endfor
01369 %!testif HAVE_ARPACK
01370 %! [v1,d1] = eigs(A, k, 'sm');
01371 %! d1 = diag(d1);
01372 %! for i=1:k
01373 %!  assert(max(abs((A - d1(i)*eye(n))*v1(:,i))),0.,1e-11)
01374 %! endfor
01375 %!testif HAVE_ARPACK
01376 %! [v1,d1] = eigs(A, k, 'lr');
01377 %! d1 = diag(d1);
01378 %! for i=1:k
01379 %!  assert(max(abs((A - d1(i)*eye(n))*v1(:,i))),0.,1e-11)
01380 %! endfor
01381 %!testif HAVE_ARPACK
01382 %! [v1,d1] = eigs(A, k, 'sr');
01383 %! d1 = diag(d1);
01384 %! for i=1:k
01385 %!  assert(max(abs((A - d1(i)*eye(n))*v1(:,i))),0.,1e-11)
01386 %! endfor
01387 %!testif HAVE_ARPACK
01388 %! [v1,d1] = eigs(A, k, 'li');
01389 %! d1 = diag(d1);
01390 %! for i=1:k
01391 %!  assert(max(abs((A - d1(i)*eye(n))*v1(:,i))),0.,1e-11)
01392 %! endfor
01393 %!testif HAVE_ARPACK
01394 %! [v1,d1] = eigs(A, k, 'si');
01395 %! d1 = diag(d1);
01396 %! for i=1:k
01397 %!  assert(max(abs((A - d1(i)*eye(n))*v1(:,i))),0.,1e-11)
01398 %! endfor
01399 
01400 */
01401 
01402 /*
01403 
01404 %% Complex hermitian tests
01405 %!shared n, k, A, d0
01406 %! n = 20;
01407 %! k = 4;
01408 %! A = full(sparse([3:n,1:n,1:(n-2)],[1:(n-2),1:n,3:n],[1i*ones(1,n-2),4*ones(1,n),-1i*ones(1,n-2)]));
01409 %! d0 = eig (A);
01410 %! [~, idx] = sort (abs(d0));
01411 %! d0 = d0(idx);
01412 %! rand("state", 42); % initialize generator to make eigs behavior reproducible
01413 %!testif HAVE_ARPACK
01414 %! d1 = eigs (A, k);
01415 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01416 %!testif HAVE_ARPACK
01417 %! d1 = eigs (A,k+1);
01418 %! assert (abs(d1), abs(d0(end:-1:(end-k))),1e-11);
01419 %!testif HAVE_ARPACK
01420 %! d1 = eigs (A, k, 'lm');
01421 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01422 %!testif HAVE_ARPACK
01423 %! d1 = eigs (A, k, 'sm');
01424 %! assert (abs(d1), abs(d0(1:k)), 1e-11);
01425 %!testif HAVE_ARPACK
01426 %! d1 = eigs (A, k, 'lr');
01427 %! [~, idx] = sort (real(abs(d0)));
01428 %! d2 = d0(idx);
01429 %! assert (real(d1), real(d2(end:-1:(end-k+1))), 1e-11);
01430 %!testif HAVE_ARPACK
01431 %! d1 = eigs (A, k, 'sr');
01432 %! [~, idx] = sort (real(abs(d0)));
01433 %! d2 = d0(idx);
01434 %! assert (real(d1), real(d2(1:k)), 1e-11);
01435 %!testif HAVE_ARPACK
01436 %! d1 = eigs (A, k, 'li');
01437 %! [~, idx] = sort (imag(abs(d0)));
01438 %! d2 = d0(idx);
01439 %! assert (sort(imag(d1)), sort(imag(d2(end:-1:(end-k+1)))), 1e-11);
01440 %!testif HAVE_ARPACK
01441 %! d1 = eigs (A, k, 'si');
01442 %! [~, idx] = sort (imag(abs(d0)));
01443 %! d2 = d0(idx);
01444 %! assert (sort(imag(d1)), sort(imag(d2(1:k))), 1e-11);
01445 %!testif HAVE_ARPACK
01446 %! d1 = eigs (A, k, 4.1);
01447 %! [~,idx0] = sort (abs(d0 - 4.1));
01448 %! [~,idx1] = sort (abs(d1 - 4.1));
01449 %! assert (abs(d1(idx1)), abs(d0(idx0(1:k))), 1e-11);
01450 %! assert (sort(imag(d1(idx1))), sort(imag(d0(idx0(1:k)))), 1e-11);
01451 %!testif HAVE_ARPACK
01452 %! d1 = eigs(A, eye(n), k, 'lm');
01453 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01454 %!testif HAVE_ARPACK
01455 %! opts.cholB=true;
01456 %! d1 = eigs(A, eye(n), k, 'lm', opts);
01457 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01458 %!testif HAVE_ARPACK
01459 %! opts.cholB=true;
01460 %! q = [2:n,1];
01461 %! opts.permB=q;
01462 %! d1 = eigs(A, eye(n)(q,q), k, 'lm', opts);
01463 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01464 %!testif HAVE_ARPACK
01465 %! opts.cholB=true;
01466 %! d1 = eigs(A, eye(n), k, 4.1, opts);
01467 %! assert (abs(abs(d1)), abs(eigs(A,k,4.1)), 1e-11);
01468 %! assert (sort(imag(abs(d1))), sort(imag(eigs(A,k,4.1))), 1e-11);
01469 %!testif HAVE_ARPACK
01470 %! opts.cholB=true;
01471 %! q = [2:n,1];
01472 %! opts.permB=q;
01473 %! d1 = eigs(A, eye(n)(q,q), k, 4.1, opts);
01474 %! assert (abs(abs(d1)), abs(eigs(A,k,4.1)), 1e-11);
01475 %! assert (sort(imag(abs(d1))), sort(imag(eigs(A,k,4.1))), 1e-11);
01476 %!testif HAVE_ARPACK
01477 %! assert (abs(eigs(A,k,4.1)), abs(eigs(A,eye(n),k,4.1)), 1e-11);
01478 %!testif HAVE_ARPACK
01479 %! assert (sort(imag(eigs(A,k,4.1))), sort(imag(eigs(A,eye(n),k,4.1))), 1e-11);
01480 %!testif HAVE_ARPACK
01481 %! fn = @(x) A * x;
01482 %! opts.issym = 0; opts.isreal = 0;
01483 %! d1 = eigs (fn, n, k, 'lm', opts);
01484 %! assert (abs(d1), abs(d0(end:-1:(end-k+1))), 1e-11);
01485 %!testif HAVE_ARPACK
01486 %! fn = @(x) A \ x;
01487 %! opts.issym = 0; opts.isreal = 0;
01488 %! d1 = eigs (fn, n, k, 'sm', opts);
01489 %! assert (abs(d1), d0(1:k), 1e-11);
01490 %!testif HAVE_ARPACK
01491 %! fn = @(x) (A - 4.1 * eye(n)) \ x;
01492 %! opts.issym = 0; opts.isreal = 0;
01493 %! d1 = eigs (fn, n, k, 4.1, opts);
01494 %! assert (abs(d1), eigs(A,k,4.1), 1e-11);
01495 %!testif HAVE_ARPACK
01496 %! [v1,d1] = eigs(A, k, 'lm');
01497 %! d1 = diag(d1);
01498 %! for i=1:k
01499 %!  assert(max(abs((A - d1(i)*eye(n))*v1(:,i))),0.,1e-11)
01500 %! endfor
01501 %!testif HAVE_ARPACK
01502 %! [v1,d1] = eigs(A, k, 'sm');
01503 %! d1 = diag(d1);
01504 %! for i=1:k
01505 %!  assert(max(abs((A - d1(i)*eye(n))*v1(:,i))),0.,1e-11)
01506 %! endfor
01507 %!testif HAVE_ARPACK
01508 %! [v1,d1] = eigs(A, k, 'lr');
01509 %! d1 = diag(d1);
01510 %! for i=1:k
01511 %!  assert(max(abs((A - d1(i)*eye(n))*v1(:,i))),0.,1e-11)
01512 %! endfor
01513 %!testif HAVE_ARPACK
01514 %! [v1,d1] = eigs(A, k, 'sr');
01515 %! d1 = diag(d1);
01516 %! for i=1:k
01517 %!  assert(max(abs((A - d1(i)*eye(n))*v1(:,i))),0.,1e-11)
01518 %! endfor
01519 %!testif HAVE_ARPACK
01520 %! [v1,d1] = eigs(A, k, 'li');
01521 %! d1 = diag(d1);
01522 %! for i=1:k
01523 %!  assert(max(abs((A - d1(i)*eye(n))*v1(:,i))),0.,1e-11)
01524 %! endfor
01525 %!testif HAVE_ARPACK
01526 %! [v1,d1] = eigs(A, k, 'si');
01527 %! d1 = diag(d1);
01528 %! for i=1:k
01529 %!  assert(max(abs((A - d1(i)*eye(n))*v1(:,i))),0.,1e-11)
01530 %! endfor
01531 
01532 */
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Defines