hide long namesshow long names
hide short namesshow short names
Integer type:  int32  int64  nag_int  show int32  show int32  show int64  show int64  show nag_int  show nag_int

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

NAG Toolbox: nag_eigen_real_gen_sparse_arnoldi (f02ek)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_eigen_real_gen_sparse_arnoldi (f02ek) computes selected eigenvalues and eigenvectors of a real sparse general matrix.

Syntax

[a, nconv, w, v, resid, user, ifail] = f02ek(n, a, icolzp, irowix, nev, ncv, sigma, monit, option, 'nz', nz, 'user', user)
[a, nconv, w, v, resid, user, ifail] = nag_eigen_real_gen_sparse_arnoldi(n, a, icolzp, irowix, nev, ncv, sigma, monit, option, 'nz', nz, 'user', user)

Description

nag_eigen_real_gen_sparse_arnoldi (f02ek) computes selected eigenvalues and the corresponding right eigenvectors of a real sparse general matrix A:
Awi = λi wi .  
A specified number, nev, of eigenvalues λi, or the shifted inverses μi=1/λi-σ, may be selected either by largest or smallest modulus, largest or smallest real part, or, largest or smallest imaginary part. Convergence is generally faster when selecting larger eigenvalues, smaller eigenvalues can always be selected by choosing a zero inverse shift (σ=0.0). When eigenvalues closest to a given real value are required then the shifted inverses of largest magnitude should be selected with shift equal to the required real value.
Note that even though A is real, λi and wi may be complex. If wi is an eigenvector corresponding to a complex eigenvalue λi, then the complex conjugate vector w-i is the eigenvector corresponding to the complex conjugate eigenvalue λ-i. The eigenvalues in a complex conjugate pair λi and λ-i are either both selected or both not selected.
The sparse matrix A is stored in compressed column storage (CCS) format. See Compressed column storage (CCS) format in the F11 Chapter Introduction.
nag_eigen_real_gen_sparse_arnoldi (f02ek) uses an implicitly restarted Arnoldi iterative method to converge approximations to a set of required eigenvalues and corresponding eigenvectors. Further algorithmic information is given in Further Comments while a fuller discussion is provided in the F12 Chapter Introduction. If shifts are to be performed then operations using shifted inverse matrices are performed using a direct sparse solver; further information on the solver used is provided in the F11 Chapter Introduction.

References

Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore
Lehoucq R B, Sorensen D C and Yang C (1998) ARPACK Users' Guide: Solution of Large-scale Eigenvalue Problems with Implicitly Restarted Arnoldi Methods SIAM, Philidelphia

Parameters

Compulsory Input Parameters

1:     n int64int32nag_int scalar
n, the order of the matrix A.
Constraint: n0.
2:     anz – double array
The array of nonzero elements (and diagonal elements if a nonzero inverse shift is to be applied) of the n by n general matrix A.
3:     icolzpn+1 int64int32nag_int array
icolzpi contains the index in a of the start of column i, for i=1,2,,n; icolzpn+1 must contain the value nz+1. Thus the number of nonzero elements in column i of A is icolzpi+1-icolzpi; when shifts are applied this includes diagonal elements irrespective of value. See Compressed column storage (CCS) format in the F11 Chapter Introduction.
4:     irowixnz int64int32nag_int array
irowixi contains the row index for each entry in a. See Compressed column storage (CCS) format in the F11 Chapter Introduction.
5:     nev int64int32nag_int scalar
The number of eigenvalues to be computed.
Constraint: 0<nev<n-1.
6:     ncv int64int32nag_int scalar
The dimension of the array w. the number of Arnoldi basis vectors to use during the computation.
At present there is no a priori analysis to guide the selection of ncv relative to nev. However, it is recommended that ncv2×nev+1. If many problems of the same type are to be solved, you should experiment with increasing ncv while keeping nev fixed for a given test problem. This will usually decrease the required number of matrix-vector operations but it also increases the work and storage required to maintain the orthogonal basis vectors. The optimal ‘cross-over’ with respect to CPU time is problem dependent and must be determined empirically.
Constraint: nev+1<ncvn.
7:     sigma – double scalar
If the Shifted Inverse Real mode has been selected then sigma contains the real shift used; otherwise sigma is not referenced. This mode can be selected by setting the appropriate options in the user-supplied function option.
8:     monit – function handle or string containing name of m-file
monit is used to monitor the progress of nag_eigen_real_gen_sparse_arnoldi (f02ek). monit may be the dummy function nag_eigen_arnoldi_monit_gen (f02ekz) if no monitoring is actually required. (nag_eigen_arnoldi_monit_gen (f02ekz) is included in the NAG Toolbox.) monit is called after the solution of each eigenvalue sub-problem and also just prior to return from nag_eigen_real_gen_sparse_arnoldi (f02ek).
[istat, user] = monit(ncv, niter, nconv, w, rzest, istat, user)

Input Parameters

1:     ncv int64int32nag_int scalar
The dimension of the arrays w and rzest. The number of Arnoldi basis vectors used during the computation.
2:     niter int64int32nag_int scalar
The number of the current Arnoldi iteration.
3:     nconv int64int32nag_int scalar
The number of converged eigenvalues so far.
4:     wncv – complex array
The first nconv elements of w contain the converged approximate eigenvalues.
5:     rzestncv – double array
The first nconv elements of rzest contain the Ritz estimates (error bounds) on the converged approximate eigenvalues.
6:     istat int64int32nag_int scalar
Set to zero.
7:     user – Any MATLAB object
monit is called from nag_eigen_real_gen_sparse_arnoldi (f02ek) with the object supplied to nag_eigen_real_gen_sparse_arnoldi (f02ek).

Output Parameters

1:     istat int64int32nag_int scalar
If set to a nonzero value nag_eigen_real_gen_sparse_arnoldi (f02ek) returns immediately with ifail=9.
2:     user – Any MATLAB object
9:     option – function handle or string containing name of m-file
You can supply non-default options to the Arnoldi eigensolver by repeated calls to nag_sparseig_real_option (f12ad) from within option. (Please note that it is only necessary to call nag_sparseig_real_option (f12ad); no call to nag_sparseig_real_init (f12aa) is required from within option.) For example, you can set the mode to Shifted Inverse Real, you can increase the Iteration Limit beyond its default and you can print varying levels of detail on the iterative process using Print Level.
If only the default options (including that the eigenvalues of largest magnitude are sought) are to be used then option may be the dummy function nag_eigen_arnoldi_option (f02eky) (nag_eigen_arnoldi_option (f02eky) is included in the NAG Toolbox). See Example for an example of using option to set some non-default options.
[icomm, comm, istat, user] = option(icomm, comm, istat, user)

Input Parameters

1:     icomm: int64int32nag_int array
Contains details of the default option set. This array must be passed as argument icomm in any call to nag_sparseig_real_option (f12ad).
2:     comm: – double array
Contains details of the default option set. This array must be passed as argument comm in any call to nag_sparseig_real_option (f12ad).
3:     istat int64int32nag_int scalar
Set to zero.
4:     user – Any MATLAB object
option is called from nag_eigen_real_gen_sparse_arnoldi (f02ek) with the object supplied to nag_eigen_real_gen_sparse_arnoldi (f02ek).

Output Parameters

1:     icomm: int64int32nag_int array
Contains data on the current options set which may be altered from the default set via calls to nag_sparseig_real_option (f12ad).
2:     comm: – double array
Contains data on the current options set which may be altered from the default set via calls to nag_sparseig_real_option (f12ad).
3:     istat int64int32nag_int scalar
If set to a nonzero value nag_eigen_real_gen_sparse_arnoldi (f02ek) returns immediately with ifail=10.
4:     user – Any MATLAB object

Optional Input Parameters

1:     nz int64int32nag_int scalar
Default: the dimension of the array a and the dimension of the array irowix. (An error is raised if these dimensions are not equal.)
The dimension of the array a and the number of nonzero elements of the matrix A and, if a nonzero shifted inverse is to be applied, all diagonal elements. Each nonzero is counted once in the latter case.
Constraint: 0nzn2.
2:     user – Any MATLAB object
user is not used by nag_eigen_real_gen_sparse_arnoldi (f02ek), but is passed to monit and option. Note that for large objects it may be more efficient to use a global variable which is accessible from the m-files than to use user.

Output Parameters

1:     anz – double array
If a nonzero shifted inverse is to be applied then the diagonal elements of A have the shift value, as supplied in sigma, subtracted.
2:     nconv int64int32nag_int scalar
The number of converged approximations to the selected eigenvalues. On successful exit, this will normally be either nev or nev+1 depending on the number of complex conjugate pairs of eigenvalues returned.
3:     wncv – complex array
The first nconv elements contain the converged approximations to the selected eigenvalues. Since complex conjugate pairs of eigenvalues appear together, it is possible (given an odd number of converged real eigenvalues) for nag_eigen_real_gen_sparse_arnoldi (f02ek) to return one more eigenvalue than requested.
4:     vldv: – double array
The first dimension of the array v will be n.
The second dimension of the array v will be ncv.
Contains the eigenvectors associated with the eigenvalue λi, for i=1,2,,nconv (stored in w). For a real eigenvalue, λj, the corresponding eigenvector is real and is stored in vij, for i=1,2,,n. For complex conjugate pairs of eigenvalues, wj+1=wj-, the real and imaginary parts of the corresponding eigenvectors are stored, respectively, in vij and vij, for i=1,2,,n. The imaginary parts stored are for the first of the conjugate pair of eigenvectors; the other eigenvector in the pair is obtained by negating these imaginary parts.
5:     residnev+1 – double array
The residual A wi - λi wi 2 for the estimates to the eigenpair λi and wi is returned in residi, for i=1,2,,nconv.
6:     user – Any MATLAB object
7:     ifail int64int32nag_int scalar
ifail=0 unless the function detects an error (see Error Indicators and Warnings).

Error Indicators and Warnings

Errors or warnings detected by the function:
   ifail=1
Constraint: n>0.
   ifail=2
Constraint: nz>0.
Constraint: nzn×n.
   ifail=3
On entry, in shifted inverse mode, the jth diagonal element of A is not defined.
   ifail=4
Constraint: icolzp1=1.
Constraint: icolzpi<icolzpi+1.
Constraint: icolzpn+1=nz+1.
   ifail=5
Constraint: irowixj<irowixj+1.
   ifail=6
Constraint: nev>0.
   ifail=7
Constraint: ncv>nev+1.
Constraint: ncvn.
   ifail=8
On entry, the matrix A-σ×I is nearly numerically singular and could not be inverted.
On entry, the matrix A-σ×I is numerically singular and could not be inverted. Try perturbing the value of σ.
   ifail=9
User requested termination in monit.
   ifail=10
User requested termination in option.
   ifail=14
Constraint: ldvn.
   ifail=21
The maximum number of iterations 0, the optional parameter Iteration Limit has been set.
   ifail=22
An internal call to nag_sparseig_real_iter (f12ab) returned with ifail=2.
This error should not occur. Please contact NAG.
   ifail=23
An internal call to nag_sparseig_real_iter (f12ab) returned with ifail=3.
   ifail=24
The maximum number of iterations has been reached.
   ifail=25
No shifts could be applied during a cycle of the implicitly restarted Arnoldi iteration.
   ifail=26
Could not build an Arnoldi factorization.
   ifail=27
Error in internal call to compute eigenvalues and corresponding error bounds of the current upper Hessenberg matrix.
Please contact NAG.
   ifail=32
An internal call to nag_sparseig_real_proc (f12ac) returned with ifail=2.
   ifail=33
The number of eigenvalues found to sufficient accuracy is zero.
   ifail=34
Internal inconsistency in the number of converged Ritz values.
   ifail=35
During calculation of a real Schur form, there was a failure to compute _ eigenvalues in a total of _ iterations.
   ifail=36
The computed Schur form could not be reordered by an internal call.
   ifail=37
In calculating eigenvectors, an internal call returned with an error.
Please contact NAG.
   ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
   ifail=-399
Your licence key may have expired or may not have been installed correctly.
   ifail=-999
Dynamic memory allocation failed.

Accuracy

The relative accuracy of a Ritz value (eigenvalue approximation), λ, is considered acceptable if its Ritz estimate Tolerance×λ. The default value for Tolerance is the machine precision given by nag_machine_precision (x02aj). The Ritz estimates are available via the monit function at each iteration in the Arnoldi process, or can be printed by setting option Print Level to a positive value.

Further Comments

nag_eigen_real_gen_sparse_arnoldi (f02ek) calls functions based on the ARPACK suite in Chapter F12. These functions use an implicitly restarted Arnoldi iterative method to converge to approximations to a set of required eigenvalues (see the F12 Chapter Introduction).
In the default Regular mode, only matrix-vector multiplications are performed using the sparse matrix A during the Arnoldi process. Each iteration is therefore cheap computationally, relative to the alternative, Shifted Inverse Real, mode described below. It is most efficient (i.e., the total number of iterations required is small) when the eigenvalues of largest magnitude are sought and these are distinct.
Although there is an option for returning the smallest eigenvalues using this mode (see Smallest Magnitude option), the number of iterations required for convergence will be far greater or the method may not converge at all. However, where convergence is achieved, Regular mode may still prove to be the most efficient since no inversions are required. Where smallest eigenvalues are sought and Regular mode is not suitable, or eigenvalues close to a given real value are sought, the Shifted Inverse Real mode should be used.
If the Shifted Inverse Real mode is used (via a call to nag_sparseig_real_option (f12ad) in option) then the matrix A-σI is used in linear system solves by the Arnoldi process. This is first factorized internally using the direct LU factorization function nag_sparse_direct_real_gen_lu (f11me). The condition number of A-σI is then calculated by a call to nag_sparse_direct_real_gen_cond (f11mg). If the condition number is too big then the matrix is considered to be nearly singular, i.e., σ is an approximate eigenvalue of A, and the function exits with an error. In this situation it is normally sufficient to perturb σ by a small amount and call nag_eigen_real_gen_sparse_arnoldi (f02ek) again. After successful factorization, subsequent solves are performed by calls to nag_sparse_direct_real_gen_solve (f11mf).
Finally, nag_eigen_real_gen_sparse_arnoldi (f02ek) transforms the eigenvectors. Each eigenvector w (real or complex) is normalized so that w2=1, and the element of largest absolute value is real.
The monitoring function monit provides some basic information on the convergence of the Arnoldi iterations. Much greater levels of detail on the Arnoldi process are available via option Print Level. If this is set to a positive value then information will be printed, by default, to standard output. The Monitoring option may be used to select a monitoring file by setting the option to a file identification (unit) number associated with Monitoring (see nag_file_open (x04ac)).

Example

This example computes the four eigenvalues of the matrix A which lie closest to the value σ=5.5 on the real line, and their corresponding eigenvectors, where A is the tridiagonal matrix with elements
aij = 2+i, j=i 3, j = i-1 -1 + ρ / 2n+2, j = i+1 ​ with ​ ρ = 10.0.  
function f02ek_example


fprintf('f02ek example results\n\n');


% This example demonstrates the use of f02ek to evaluate a number of
% eigenvalues of a sparse matrix closest to a given value.

% Use existing sparse matrix.
load('west0479.mat')
W = west0479;

[irow,icol,a] = find(W);
n = int64(size(W,1));
nnz = int64(size(irow,1));

% Add zero diagonals to be summed
irow(nnz+1:nnz+n) = (1:n);
icol(nnz+1:nnz+n) = (1:n);
a(nnz+1:nnz+n) = 0.0;
nnz = nnz + n;

% Use f11za to get Compressed column storage array icolzp.
irow = int64(irow);
icol = int64(icol);
dup = 'S';
zero = 'K';
[nnz, a, icol, irow, icolzp, ifail] = ...
  f11za(...
        n, nnz, a, icol, irow, dup, zero);

% Evaluate nev eigenvalues (w) closest to sigma.
nev = int64(20);
ncv = int64(60);
sigma = 4.0;
% Set some options via user array and function argument option.
% user(1) = print level, iuser(2) = iteration limit,
% user(3)>0 means shifted-invert mode
% user(4)>0 means print monitoring info
user = [0, 500, 1, 0];

[a, nconv, w, v, resid, user, ifail] = ...
  f02ek(...
        n, a, icolzp, irow, nev, ncv, sigma, @monit, @option, 'user', user);

fprintf('\n The %d Ritz values closest to %8.2e are:\n', nconv, sigma);
disp(w(1:nconv));

fig1 = figure;
plot(w(1:nconv),'k+');
title('Eigenvalues of West0479 closest to (4,0)');
xlabel('real');
ylabel('imag');


function [istat, user] = monit(ncv, niter, nconv, w, rzest, istat, user);
  if (user(4)>0)
    if (niter==1 && user(3)>0)
      fprintf('The following Ritz values (mu) are related to the\n');
      fprintf('true eigenvalues (lambda) by lambda = sigma + 1/mu\n\n');
    end
    fprintf('Iteration number %d\n', niter);
    fprintf(' Ritz values converged so far (%d):\n', nconv);
    for i = 1:nconv
      fprintf(' %4d (%13.5e,%13.5e)\n',  i, real(w(i)), imag(w(i)));
    end
    fprintf(' Next (unconverged) Ritz value:\n');
    fprintf(' %4d (%13.5e,%13.5e)\n', nconv+1, real(w(nconv+1)), ...
            imag(w(nconv+1)));
  end
  istat = int64(0);

function [icomm, comm, istat, user] = option(icomm, comm, istat, user);
  if (user(1)>0)
    [icomm, comm, ifail] = ...
      f12ad(strcat('Print Level=', int2str(user(1))), icomm, comm);
    istat = max(istat,ifail);
  end
  if (user(2)>100)
    [icomm, comm, ifail] = ...
      f12ad(strcat('Iteration Limit=', int2str(user(2))), icomm, comm);
    istat = max(istat,ifail);
  end
  if (user(3)>0)
    [icomm, comm, ifail] = ...
      f12ad('Shifted Inverse Real', icomm, comm);
    istat = max(istat,ifail);
  end
f02ek example results


 The 20 Ritz values closest to 4.00e+00 are:
   3.9852 + 0.0000i
   4.7459 + 0.0000i
   3.1237 - 0.1509i
   3.1237 + 0.1509i
   4.3290 - 1.0151i
   4.3290 + 1.0151i
   2.8651 - 0.6759i
   2.8651 + 0.6759i
   2.5930 + 0.0000i
   5.4060 + 0.0000i
   5.8229 + 0.0000i
   2.3494 - 0.8839i
   2.3494 + 0.8839i
   2.2071 - 0.3953i
   2.2071 + 0.3953i
   2.5364 - 1.4458i
   2.5364 + 1.4458i
   2.0675 + 0.0000i
   1.9248 - 0.2785i
   1.9248 + 0.2785i

f02ek_fig1.png

Optional Parameters

Internally nag_eigen_real_gen_sparse_arnoldi (f02ek) calls functions from the suite nag_sparseig_real_init (f12aa), nag_sparseig_real_iter (f12ab), nag_sparseig_real_proc (f12ac), nag_sparseig_real_option (f12ad) and nag_sparseig_real_monit (f12ae). Several optional parameters for these computational functions define choices in the problem specification or the algorithm logic. In order to reduce the number of formal arguments of nag_eigen_real_gen_sparse_arnoldi (f02ek) these optional parameters are also used here and have associated default values that are usually appropriate. Therefore, you need only specify those optional parameters whose values are to be different from their default values.
Optional parameters may be specified via the user-supplied function option in the call to nag_eigen_real_gen_sparse_arnoldi (f02ek). option must be coded such that one call to nag_sparseig_real_option (f12ad) is necessary to set each optional parameter. All optional parameters you do not specify are set to their default values.
The remainder of this section can be skipped if you wish to use the default values for all optional parameters.
The following is a list of the optional parameters available. A full description of each optional parameter is provided in Description of the s.

Description of the Optional Parameters

For each option, we give a summary line, a description of the optional parameter and details of constraints.
The summary line contains:
Keywords and character values are case and white space insensitive.
Advisory  i
Default =0  
If the optional parameter List is set then optional parameter specifications are listed in a List file by setting the option to a file identification (unit) number associated with Advisory messages (see nag_file_set_unit_advisory (x04ab) and nag_file_open (x04ac)).
Defaults  
This special keyword may be used to reset all optional parameters to their default values.
Iteration Limit  i
Default = 300  
The limit on the number of Arnoldi iterations that can be performed before nag_eigen_real_gen_sparse_arnoldi (f02ek) exits with ifail0.
Largest Magnitude  
Default
Largest Imaginary  
Largest Real  
Smallest Imaginary  
Smallest Magnitude  
Smallest Real  
The Arnoldi iterative method converges on a number of eigenvalues with given properties. The default is to compute the eigenvalues of largest magnitude using Largest Magnitude. Alternatively, eigenvalues may be chosen which have Largest Real part, Largest Imaginary part, Smallest Magnitude, Smallest Real part or Smallest Imaginary part.
Note that these options select the eigenvalue properties for eigenvalues of OP the linear operator determined by the computational mode and problem type.
Nolist  
Default
List  
Normally each optional parameter specification is not printed to the advisory channel as it is supplied. Optional parameter List may be used to enable printing and optional parameter Nolist may be used to suppress the printing.
Monitoring  i
Default = -1
If i>0, monitoring information is output to channel number i during the solution of each problem; this may be the same as the Advisory channel number. The type of information produced is dependent on the value of Print Level, see the description of the optional parameter Print Level for details of the information produced. Please see nag_file_open (x04ac) to associate a file with a given channel number.
Print Level  i
Default = 0
This controls the amount of printing produced by nag_eigen_real_gen_sparse_arnoldi (f02ek) as follows.
=0 No output except error messages.
>0 The set of selected options.
=2 Problem and timing statistics when all calls to nag_sparseig_real_iter (f12ab) have been completed.
5 A single line of summary output at each Arnoldi iteration.
10 If Monitoring>0, then at each iteration, the length and additional steps of the current Arnoldi factorization and the number of converged Ritz values; during re-orthogonalization, the norm of initial/restarted starting vector.
20 Problem and timing statistics on final exit from nag_sparseig_real_iter (f12ab). If Monitoring>0, then at each iteration, the number of shifts being applied, the eigenvalues and estimates of the Hessenberg matrix H, the size of the Arnoldi basis, the wanted Ritz values and associated Ritz estimates and the shifts applied; vector norms prior to and following re-orthogonalization.
30 If Monitoring>0, then on final iteration, the norm of the residual; when computing the Schur form, the eigenvalues and Ritz estimates both before and after sorting; for each iteration, the norm of residual for compressed factorization and the compressed upper Hessenberg matrix H; during re-orthogonalization, the initial/restarted starting vector; during the Arnoldi iteration loop, a restart is flagged and the number of the residual requiring iterative refinement; while applying shifts, the indices of the shifts being applied.
40 If Monitoring>0, then during the Arnoldi iteration loop, the Arnoldi vector number and norm of the current residual; while applying shifts, key measures of progress and the order of H; while computing eigenvalues of H, the last rows of the Schur and eigenvector matrices; when computing implicit shifts, the eigenvalues and Ritz estimates of H.
50 If Monitoring>0, then during Arnoldi iteration loop: norms of key components and the active column of H, norms of residuals during iterative refinement, the final upper Hessenberg matrix H; while applying shifts: number of shifts, shift values, block indices, updated matrix H; while computing eigenvalues of H: the matrix H, the computed eigenvalues and Ritz estimates.
Regular  
Default
Shifted Inverse Real  
These options define the computational mode which in turn defines the form of operation OPx to be performed.
Given a standard eigenvalue problem in the form Ax=λx then the following modes are available with the appropriate operator OPx.
Regular OP=A
Shifted Inverse Real OP=A-σI-1 where σ is real
Tolerance  r
Default = ε  
An approximate eigenvalue has deemed to have converged when the corresponding Ritz estimate is within Tolerance relative to the magnitude of the eigenvalue.

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015