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_sparseig_real_symm_monit (f12fe)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_sparseig_real_symm_monit (f12fe) can be used to return additional monitoring information during computation. It is in a suite of functions which includes nag_sparseig_real_symm_init (f12fa), nag_sparseig_real_symm_iter (f12fb), nag_sparseig_real_symm_proc (f12fc) and nag_sparseig_real_symm_option (f12fd).

Syntax

[niter, nconv, ritz, rzest] = f12fe(icomm, comm)
[niter, nconv, ritz, rzest] = nag_sparseig_real_symm_monit(icomm, comm)

Description

The suite of functions is designed to calculate some of the eigenvalues, λ , (and optionally the corresponding eigenvectors, x ) of a standard eigenvalue problem Ax = λx , or of a generalized eigenvalue problem Ax = λBx  of order n , where n  is large and the coefficient matrices A  and B  are sparse, real and symmetric. The suite can also be used to find selected eigenvalues/eigenvectors of smaller scale dense, real and symmetric problems.
On an intermediate exit from nag_sparseig_real_symm_iter (f12fb) with irevcm = 4 , nag_sparseig_real_symm_monit (f12fe) may be called to return monitoring information on the progress of the Arnoldi iterative process. The information returned by nag_sparseig_real_symm_monit (f12fe) is:
the number of the current Arnoldi iteration;
the number of converged eigenvalues at this point;
the real and imaginary parts of the converged eigenvalues;
the error bounds on the converged eigenvalues.
nag_sparseig_real_symm_monit (f12fe) does not have an equivalent function from the ARPACK package which prints various levels of detail of monitoring information through an output channel controlled via a argument value (see Lehoucq et al. (1998) for details of ARPACK routines). nag_sparseig_real_symm_monit (f12fe) should not be called at any time other than immediately following an irevcm = 4  return from nag_sparseig_real_symm_iter (f12fb).

References

Lehoucq R B (2001) Implicitly restarted Arnoldi methods and subspace iteration SIAM Journal on Matrix Analysis and Applications 23 551–562
Lehoucq R B and Scott J A (1996) An evaluation of software for computing eigenvalues of sparse nonsymmetric matrices Preprint MCS-P547-1195 Argonne National Laboratory
Lehoucq R B and Sorensen D C (1996) Deflation techniques for an implicitly restarted Arnoldi iteration SIAM Journal on Matrix Analysis and Applications 17 789–821
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:     icomm: int64int32nag_int array
The dimension of the array icomm must be at least max1,licomm, where licomm is passed to the setup function  (see nag_sparseig_real_symm_init (f12fa))
The array icomm output by the preceding call to nag_sparseig_real_symm_iter (f12fb).
2:     comm: – double array
The dimension of the array comm must be at least max1, 3×n+ ncv× ncv+ 8×ncv +60  (see nag_sparseig_real_symm_init (f12fa))
The array comm output by the preceding call to nag_sparseig_real_symm_iter (f12fb).

Optional Input Parameters

None.

Output Parameters

1:     niter int64int32nag_int scalar
The number of the current Arnoldi iteration.
2:     nconv int64int32nag_int scalar
The number of converged eigenvalues so far.
3:     ritz: – double array
The dimension of the array ritz will be ncv (see nag_sparseig_real_symm_init (f12fa))
The first nconv locations of the array ritz contain the real converged approximate eigenvalues.
4:     rzest: – double array
The dimension of the array rzest will be ncv (see nag_sparseig_real_symm_init (f12fa))
The first nconv locations of the array rzest contain the Ritz estimates (error bounds) on the real nconv converged approximate eigenvalues.

Error Indicators and Warnings

None.

Accuracy

A Ritz value, λ , is deemed to have converged if its Ritz estimate Tolerance × λ . The default Tolerance used is the machine precision given by nag_machine_precision (x02aj).

Further Comments

None.

Example

This example solves Kx = λKGx  using the Buckling option (see nag_sparseig_real_symm_option (f12fd), where K  and KG  are obtained by the finite element method applied to the one-dimensional discrete Laplacian operator 2u x2  on 0,1 , with zero Dirichlet boundary conditions using piecewise linear elements. The shift, σ , is a real number, and the operator used in the Buckling iterative process is OP = invK - σKG × K  and B = K .
function f12fe_example


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

n = int64(100);
nev = int64(4);
ncv = int64(10);

irevcm = int64(0);
resid = zeros(n,1);
v     = zeros(n,ncv);
x     = zeros(n,1);
mx    = zeros(n,1);
imon = 1;

sigma = 1;

% Setup and factorize A - sigma*B
h   = 1/double(n+1);
ad(1:n)  =  2/h - sigma*4*h/6;
adl(1:n) = -1/h - sigma*h/6;
adu(1:n) = adl(1:n);

[adl, ad, adu, adu2, ipiv, info] = f07cd( ...
                                          adl, ad, adu);

% Initialisation Step
[icomm, comm, ifail] = f12fa( ...
                              n, nev, ncv);

% Set Optional Parameters
[icomm, comm, ifail] = f12fd( ...
                              'Generalized', icomm, comm);
[icomm, comm, ifail] = f12fd( ...
                              'Buckling', icomm, comm);

% Solve
while (irevcm ~= 5)
  [irevcm, resid, v, x, mx, nshift, comm, icomm, ifail] = ...
    f12fb( ...
           irevcm, resid, v, x, mx, comm, icomm);
   if (irevcm == -1)
    % Solve (A-sigma*B)y = Ax
    mx = f12fe_Ax(n, x);
    [x, info] = f07ce( ...
                       'N', adl, ad, adu, adu2, ipiv, mx);
  elseif (irevcm == 1)
    % Solve (A-sigma*B)y = Ax, Ax in mx
     [x, info] = f07ce( ...
                       'N', adl, ad, adu, adu2, ipiv, mx);
  elseif (irevcm == 2)
    % y = Ax
    mx = f12fe_Ax(n, x);
  elseif (irevcm == 4 && imon == 1)
    [niter, nconv, ritz, rzest] = ...
    f12fe(icomm, comm);

    fprintf(['Iteration %2d, No. converged = %d, ', ...
             'norm of estimates = %10.2e\n'], ...
            niter, nconv, norm(rzest(1:nev),2));
  end
end

% Post-process to compute eigenvalues/vectors
[nconv, d, z, v, comm, icomm, ifail] = ...
f12fc( ...
       sigma, resid, v, comm, icomm);

fprintf('\nThe %d Eigenvalues closest to %7.2f are:\n',nconv, sigma);
disp(d(1:nconv));



function [y] = f12fe_Ax(n,x)

  y = zeros(n,1);
  h = 1/double(n+1);

  dd = 2;
  dl = -1;
  du = -1;

  y(1) = dd*x(1) + du*x(2);
  for j=2:n-1
    y(j) = dl*x(j-1) + dd*x(j) + du*x(j+1);
  end
  y(n) = dl*x(n-1) + dd*x(n);
  y = y/h;
f12fe example results

Iteration  1, No. converged = 0, norm of estimates =   2.05e-06
Iteration  2, No. converged = 2, norm of estimates =   6.08e-11
Iteration  3, No. converged = 3, norm of estimates =   5.27e-15

The 4 Eigenvalues closest to    1.00 are:
    9.8704
   39.4912
   88.8909
  158.1175


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