Home > tt2 > exp > lobpcg.m

lobpcg

PURPOSE ^

LOBPCG solves Hermitian partial eigenproblems using preconditioning

SYNOPSIS ^

function [blockVectorX,lambda,varargout] =lobpcg(blockVectorX,operatorA,varargin)

DESCRIPTION ^

LOBPCG solves Hermitian partial eigenproblems using preconditioning

 [blockVectorX,lambda]=lobpcg(blockVectorX,operatorA)

 outputs the array of algebraic smallest eigenvalues lambda and
 corresponding matrix of orthonormalized eigenvectors blockVectorX of the
 Hermitian (full or sparse) operator operatorA using input matrix
 blockVectorX as an initial guess, without preconditioning, somewhat
 similar to 

 opts.issym=1;opts.isreal=1;K=size(blockVectorX,2);
 [blockVectorX,lambda]=eigs(operatorA,K,'SR',opts);

 for real symmetric operator operatorA, or

 K=size(blockVectorX,2);[blockVectorX,lambda]=eigs(operatorA,K,'SR');
 for Hermitian operator operatorA. 

 [blockVectorX,lambda,failureFlag]=lobpcg(blockVectorX,operatorA) 
 also returns a convergence flag.  
 If failureFlag is 0 then all the eigenvalues converged; otherwise not all
 converged.

 [blockVectorX,lambda,failureFlag,lambdaHistory,residualNormsHistory]=...
 lobpcg(blockVectorX,'operatorA','operatorB','operatorT',blockVectorY,...
 residualTolerance,maxIterations,verbosityLevel);

 computes smallest eigenvalues lambda and corresponding eigenvectors
 blockVectorX of the generalized eigenproblem Ax=lambda Bx, where 
 Hermitian operators operatorA and operatorB are given as functions, as
 well as a preconditioner, operatorT. The operators operatorB and
 operatorT must be in addition POSITIVE DEFINITE. To compute the largest
 eigenpairs of operatorA, simply apply the code to operatorA multiplied by
 -1. The code does not involve ANY matrix factorizations of operratorA and
 operatorB, thus, e.g., it preserves the sparsity and the structure of
 operatorA and operatorB. 

 residualTolerance and maxIterations control tolerance and max number of
 steps, and verbosityLevel = 0, 1, or 2 controls the amount of printed
 info. lambdaHistory is a matrix with all iterative lambdas, and
 residualNormsHistory are matrices of the history of 2-norms of residuals

 Required input: 
   blockVectorX - initial approximation to eigenvectors, full or sparse
   matrix n-by-blockSize. blockVectorX must be full rank. 
   operatorA - the operator of the problem, can be given as a matrix or as
   an M-file 

 Optional function input:
   operatorB - the second operator, if solving a generalized eigenproblem, 
       can be given as a matrix or as an M-file; by default, or if empty,
       operatorB=I.
   operatorT - preconditioner, must be given by an M-file; by default,
   operatorT=I.

 Optional constraints input: 
   blockVectorY - a full or sparse n-by-sizeY matrix of constraints, where
   sizeY < n. The iterations will be performed in the (operatorB-)
   orthogonal complement of the column-space of blockVectorY. blockVectorY
   must be full rank. 

 Optional scalar input parameters:
   residualTolerance - tolerance, by default,
   residualTolerance=n*sqrt(eps) maxIterations - max number of iterations,
   by default, maxIterations = min(n,20) verbosityLevel - either 0 (no
   info), 1, or 2 (with pictures); by default, verbosityLevel = 0.

 Required output: blockVectorX and lambda are computed blockSize
 eigenpairs, where blockSize=size(blockVectorX,2) for the initial guess
 blockVectorX if it is full rank.  

 Optional output: failureFlag, lambdaHistory and residualNormsHistory are
 described above.

 Functions operatorA(blockVectorX), operatorB(blockVectorX) and
 operatorT(blockVectorX) must support blockVectorX being a matrix, not
 just a column vector.

 Every iteration involves one application of operatorA and operatorB, and
 one of operatorT. 

 Main memory requirements: 6 (9 if isempty(operatorB)=0) matrices of the
 same size as blockVectorX, 2 matrices of the same size as blockVectorY
 (if present), and two square matrices of the size 3*blockSize. 

 The following
 Example:

 operatorA = 100.*delsq(numgrid('S',21)); [n,n]=size(operatorA);
 [blockVectorX,lambda,failureFlag]=lobpcg(randn(n,10),operatorA,1e-5,50,2);

 attempts to compute 10 first eigenpairs of the Poisson operator 
 in a 2x2 square with the mesh size 1/10 without preconditioning,
 but not all eigenpairs converge after 50 steps, so failureFlag=1.  

 The next 
 Example:

 operatorA = 100.*delsq(numgrid('S',21)); [n,n]=size(operatorA);
 blockVectorY=[];lambda_all=[];
 for j=1:5
   [blockVectorX,lambda]=...
                    lobpcg(randn(n,2),operatorA,blockVectorY,1e-5,200,2);
   blockVectorY=[blockVectorY,blockVectorX];
   lambda_all=[lambda_all' lambda']';
 end

 attemps to compute the same eigenpairs by calling the code 5 times 
 with blockSize=2 using orthogonalization to the previously founded
 eigenvectors.  

The following M-script:

 global R_cholinc
 operatorA = 100.*delsq(numgrid('S',21)); [n,n]=size(operatorA);
 R_cholinc=cholinc(operatorA,1e-3);
 [blockVectorX,lambda,failureFlag]=...
                   lobpcg(randn(n,10),operatorA,[],'precsol',1e-5,50,2);
 
 computes the same eigenpairs in less then 25 steps, so that failureFlag=0
 using the preconditioner function precsol:

 function blockVectorX=precsol(V)
 global R_cholinc 
 blockVectorX=R_cholinc\(R_cholinc'\V);
 In this example, operatorB=[] must be present in the input parameters. 
 [blockVectorX,lambda,failureFlag]=...
              lobpcg(randn(n,10),operatorA,speye(n),'precsol',1e-5,50,2);

produces similar answers, but is somewhat slower and needs more memory.  

%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

 This main function LOBPCG is a version of 
 the preconditioned conjugate gradient method (Algorithm 5.1) described in
 A. V. Knyazev, Toward the Optimal Preconditioned Eigensolver:
 Locally Optimal Block Preconditioned Conjugate Gradient Method,
 SIAM Journal on Scientific Computing 23 (2001), no. 2, pp. 517-541. 
 http://dx.doi.org/10.1137/S1064827500366124

 Known bugs/features:

 - an excessively small requested tolerance may result in often restarts
 and instability. The code is not written to produce an eps-level
 accuracy! Use common sense.  

 - the code may be very sensitive to the number of eigenpairs computed,
 if there is a cluster of eigenvalues not completely included, cf. 

 operatorA=diag([1 1.99 2:99]);
 [blockVectorX,lambda]=lobpcg(randn(100,1),operatorA,1e-10,80,2);
 [blockVectorX,lambda]=lobpcg(randn(100,2),operatorA,1e-10,80,2);
 [blockVectorX,lambda]=lobpcg(randn(100,3),operatorA,1e-10,80,2);

%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
 The main distribution site: 
 http://math.ucdenver.edu/~aknyazev/

 A C-version of this code is a part of the 
 http://code.google.com/p/blopex/
 package and is directly available, e.g., in PETSc and HYPRE.  
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function [blockVectorX,lambda,varargout] = ...
0002     lobpcg(blockVectorX,operatorA,varargin)
0003 %LOBPCG solves Hermitian partial eigenproblems using preconditioning
0004 %
0005 % [blockVectorX,lambda]=lobpcg(blockVectorX,operatorA)
0006 %
0007 % outputs the array of algebraic smallest eigenvalues lambda and
0008 % corresponding matrix of orthonormalized eigenvectors blockVectorX of the
0009 % Hermitian (full or sparse) operator operatorA using input matrix
0010 % blockVectorX as an initial guess, without preconditioning, somewhat
0011 % similar to
0012 %
0013 % opts.issym=1;opts.isreal=1;K=size(blockVectorX,2);
0014 % [blockVectorX,lambda]=eigs(operatorA,K,'SR',opts);
0015 %
0016 % for real symmetric operator operatorA, or
0017 %
0018 % K=size(blockVectorX,2);[blockVectorX,lambda]=eigs(operatorA,K,'SR');
0019 % for Hermitian operator operatorA.
0020 %
0021 % [blockVectorX,lambda,failureFlag]=lobpcg(blockVectorX,operatorA)
0022 % also returns a convergence flag.
0023 % If failureFlag is 0 then all the eigenvalues converged; otherwise not all
0024 % converged.
0025 %
0026 % [blockVectorX,lambda,failureFlag,lambdaHistory,residualNormsHistory]=...
0027 % lobpcg(blockVectorX,'operatorA','operatorB','operatorT',blockVectorY,...
0028 % residualTolerance,maxIterations,verbosityLevel);
0029 %
0030 % computes smallest eigenvalues lambda and corresponding eigenvectors
0031 % blockVectorX of the generalized eigenproblem Ax=lambda Bx, where
0032 % Hermitian operators operatorA and operatorB are given as functions, as
0033 % well as a preconditioner, operatorT. The operators operatorB and
0034 % operatorT must be in addition POSITIVE DEFINITE. To compute the largest
0035 % eigenpairs of operatorA, simply apply the code to operatorA multiplied by
0036 % -1. The code does not involve ANY matrix factorizations of operratorA and
0037 % operatorB, thus, e.g., it preserves the sparsity and the structure of
0038 % operatorA and operatorB.
0039 %
0040 % residualTolerance and maxIterations control tolerance and max number of
0041 % steps, and verbosityLevel = 0, 1, or 2 controls the amount of printed
0042 % info. lambdaHistory is a matrix with all iterative lambdas, and
0043 % residualNormsHistory are matrices of the history of 2-norms of residuals
0044 %
0045 % Required input:
0046 %   blockVectorX - initial approximation to eigenvectors, full or sparse
0047 %   matrix n-by-blockSize. blockVectorX must be full rank.
0048 %   operatorA - the operator of the problem, can be given as a matrix or as
0049 %   an M-file
0050 %
0051 % Optional function input:
0052 %   operatorB - the second operator, if solving a generalized eigenproblem,
0053 %       can be given as a matrix or as an M-file; by default, or if empty,
0054 %       operatorB=I.
0055 %   operatorT - preconditioner, must be given by an M-file; by default,
0056 %   operatorT=I.
0057 %
0058 % Optional constraints input:
0059 %   blockVectorY - a full or sparse n-by-sizeY matrix of constraints, where
0060 %   sizeY < n. The iterations will be performed in the (operatorB-)
0061 %   orthogonal complement of the column-space of blockVectorY. blockVectorY
0062 %   must be full rank.
0063 %
0064 % Optional scalar input parameters:
0065 %   residualTolerance - tolerance, by default,
0066 %   residualTolerance=n*sqrt(eps) maxIterations - max number of iterations,
0067 %   by default, maxIterations = min(n,20) verbosityLevel - either 0 (no
0068 %   info), 1, or 2 (with pictures); by default, verbosityLevel = 0.
0069 %
0070 % Required output: blockVectorX and lambda are computed blockSize
0071 % eigenpairs, where blockSize=size(blockVectorX,2) for the initial guess
0072 % blockVectorX if it is full rank.
0073 %
0074 % Optional output: failureFlag, lambdaHistory and residualNormsHistory are
0075 % described above.
0076 %
0077 % Functions operatorA(blockVectorX), operatorB(blockVectorX) and
0078 % operatorT(blockVectorX) must support blockVectorX being a matrix, not
0079 % just a column vector.
0080 %
0081 % Every iteration involves one application of operatorA and operatorB, and
0082 % one of operatorT.
0083 %
0084 % Main memory requirements: 6 (9 if isempty(operatorB)=0) matrices of the
0085 % same size as blockVectorX, 2 matrices of the same size as blockVectorY
0086 % (if present), and two square matrices of the size 3*blockSize.
0087 %
0088 % The following
0089 % Example:
0090 %
0091 % operatorA = 100.*delsq(numgrid('S',21)); [n,n]=size(operatorA);
0092 % [blockVectorX,lambda,failureFlag]=lobpcg(randn(n,10),operatorA,1e-5,50,2);
0093 %
0094 % attempts to compute 10 first eigenpairs of the Poisson operator
0095 % in a 2x2 square with the mesh size 1/10 without preconditioning,
0096 % but not all eigenpairs converge after 50 steps, so failureFlag=1.
0097 %
0098 % The next
0099 % Example:
0100 %
0101 % operatorA = 100.*delsq(numgrid('S',21)); [n,n]=size(operatorA);
0102 % blockVectorY=[];lambda_all=[];
0103 % for j=1:5
0104 %   [blockVectorX,lambda]=...
0105 %                    lobpcg(randn(n,2),operatorA,blockVectorY,1e-5,200,2);
0106 %   blockVectorY=[blockVectorY,blockVectorX];
0107 %   lambda_all=[lambda_all' lambda']';
0108 % end
0109 %
0110 % attemps to compute the same eigenpairs by calling the code 5 times
0111 % with blockSize=2 using orthogonalization to the previously founded
0112 % eigenvectors.
0113 %
0114 %The following M-script:
0115 %
0116 % global R_cholinc
0117 % operatorA = 100.*delsq(numgrid('S',21)); [n,n]=size(operatorA);
0118 % R_cholinc=cholinc(operatorA,1e-3);
0119 % [blockVectorX,lambda,failureFlag]=...
0120 %                   lobpcg(randn(n,10),operatorA,[],'precsol',1e-5,50,2);
0121 %
0122 % computes the same eigenpairs in less then 25 steps, so that failureFlag=0
0123 % using the preconditioner function precsol:
0124 %
0125 % function blockVectorX=precsol(V)
0126 % global R_cholinc
0127 % blockVectorX=R_cholinc\(R_cholinc'\V);
0128 % In this example, operatorB=[] must be present in the input parameters.
0129 % [blockVectorX,lambda,failureFlag]=...
0130 %              lobpcg(randn(n,10),operatorA,speye(n),'precsol',1e-5,50,2);
0131 %
0132 %produces similar answers, but is somewhat slower and needs more memory.
0133 %
0134 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0135 %
0136 % This main function LOBPCG is a version of
0137 % the preconditioned conjugate gradient method (Algorithm 5.1) described in
0138 % A. V. Knyazev, Toward the Optimal Preconditioned Eigensolver:
0139 % Locally Optimal Block Preconditioned Conjugate Gradient Method,
0140 % SIAM Journal on Scientific Computing 23 (2001), no. 2, pp. 517-541.
0141 % http://dx.doi.org/10.1137/S1064827500366124
0142 %
0143 % Known bugs/features:
0144 %
0145 % - an excessively small requested tolerance may result in often restarts
0146 % and instability. The code is not written to produce an eps-level
0147 % accuracy! Use common sense.
0148 %
0149 % - the code may be very sensitive to the number of eigenpairs computed,
0150 % if there is a cluster of eigenvalues not completely included, cf.
0151 %
0152 % operatorA=diag([1 1.99 2:99]);
0153 % [blockVectorX,lambda]=lobpcg(randn(100,1),operatorA,1e-10,80,2);
0154 % [blockVectorX,lambda]=lobpcg(randn(100,2),operatorA,1e-10,80,2);
0155 % [blockVectorX,lambda]=lobpcg(randn(100,3),operatorA,1e-10,80,2);
0156 %
0157 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0158 % The main distribution site:
0159 % http://math.ucdenver.edu/~aknyazev/
0160 %
0161 % A C-version of this code is a part of the
0162 % http://code.google.com/p/blopex/
0163 % package and is directly available, e.g., in PETSc and HYPRE.
0164 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0165 
0166 %   License: GNU LGPL ver 2.1 or above
0167 %   Copyright (c) 2000-2010 A.V. Knyazev, BLOPEX team
0168 %   $Revision: 4.12 $  $Date: 14-Mar-2010
0169 %   Tested in MATLAB 6.5-7.9.0.529 and Octave 3.2.3.
0170 
0171 
0172 %Begin
0173 
0174 % Check if we are in MATLAB or Octave.  ver('MATLAB') will return an empty
0175 % struct if we are in Octave. Used below to suppress the graphics.
0176 % version=ver('MATLAB');
0177 % if all(size(version)==0)
0178 %     matlabFlag=0;
0179 % else
0180 %     matlabFlag=1;
0181 % end
0182 matlabFlag=1;
0183 
0184 
0185 % constants
0186 
0187 CONVENTIONAL_CONSTRAINTS = 1;
0188 SYMMETRIC_CONSTRAINTS = 2;
0189 
0190 %Initial settings
0191 
0192 failureFlag = 1;
0193 if nargin < 2
0194     error('BLOPEX:lobpcg:NotEnoughInputs',...
0195     strcat('There must be at least 2 input agruments: ',...
0196         'blockVectorX and operatorA'));
0197 end
0198 if nargin > 8
0199     warning('BLOPEX:lobpcg:TooManyInputs',...
0200         strcat('There must be at most 8 input agruments ',...
0201         'unless arguments are passed to a function'));
0202 end
0203 
0204 if ischar(blockVectorX)
0205     error('BLOPEX:lobpcg:FirstInputString',...
0206         'The first input argument blockVectorX cannot be a string');
0207 end
0208 [n,blockSize]=size(blockVectorX);
0209 if blockSize > n
0210     error('BLOPEX:lobpcg:FirstInputFat',...
0211     'The first input argument blockVectorX must be tall, not fat');
0212 end
0213 if n < 6
0214     error('BLOPEX:lobpcg:MatrixTooSmall',...
0215         'The code does not work for matrices of small sizes');
0216 end
0217 
0218 if ischar(operatorA)
0219     nA = size(operatorA,1);
0220     if any(size(operatorA) ~= nA)
0221         error('BLOPEX:lobpcg:MatrixNotSquare',...
0222             'operatorA must be a square matrix or a string');
0223     end
0224     if size(operatorA) ~= n
0225         error('BLOPEX:lobpcg:MatrixWrongSize',...
0226         ['The size ' int2str(size(operatorA))...
0227             ' of operatorA is not the same as ' int2str(n)...
0228             ' - the number of rows of blockVectorX']);
0229     end
0230 end
0231 
0232 count_string = 0;
0233 
0234 operatorT = [];
0235 operatorB = [];
0236 residualTolerance = [];
0237 maxIterations = [];
0238 verbosityLevel = [];
0239 blockVectorY = []; sizeY = 0;
0240 for j = 1:nargin-2
0241     if isequal(size(varargin{j}),[n,n])
0242         if isempty(operatorB)
0243             operatorB = varargin{j};
0244         else
0245             error('BLOPEX:lobpcg:TooManyMatrixInputs',...
0246         strcat('Too many matrix input arguments. ',...
0247         'Preconditioner operatorT must be an M-function'));
0248         end
0249     elseif isequal(size(varargin{j},1),n) && size(varargin{j},2) < n
0250         if isempty(blockVectorY)
0251             blockVectorY = varargin{j};
0252             sizeY=size(blockVectorY,2);
0253         else
0254             error('BLOPEX:lobpcg:WrongConstraintsFormat',...
0255             'Something wrong with blockVectorY input argument');
0256         end
0257     elseif ischar(varargin{j})
0258         if count_string == 0
0259             if isempty(operatorB)
0260                 operatorB = varargin{j};
0261                 count_string = count_string + 1;
0262             else
0263                 operatorT = varargin{j};
0264             end
0265         elseif count_string == 1
0266             operatorT = varargin{j};
0267         else
0268             error('BLOPEX:lobpcg:TooManyStringInputs',...
0269                 'Too many string input arguments');
0270         end
0271     elseif isequal(size(varargin{j}),[n,n])
0272         error('BLOPEX:lobpcg:WrongPreconditionerFormat',...
0273         'Preconditioner operatorT must be an M-function');
0274     elseif max(size(varargin{j})) == 1
0275         if isempty(residualTolerance)
0276             residualTolerance = varargin{j};
0277         elseif isempty(maxIterations)
0278             maxIterations = varargin{j};
0279         elseif isempty(verbosityLevel)
0280             verbosityLevel = varargin{j};
0281         else
0282             error('BLOPEX:lobpcg:TooManyScalarInputs',...
0283                 'Too many scalar parameters, need only three');
0284         end
0285     elseif isempty(varargin{j})
0286         if isempty(operatorB)
0287             count_string = count_string + 1;
0288         elseif ~isempty(operatorT)
0289             count_string = count_string + 1;
0290         elseif ~isempty(blockVectorY)
0291             error('BLOPEX:lobpcg:UnrecognizedEmptyInput',...
0292                ['Unrecognized empty input argument number ' int2str(j+2)]);
0293         end
0294     else
0295         error('BLOPEX:lobpcg:UnrecognizedInput',...
0296             ['Input argument number ' int2str(j+2) ' not recognized.']);
0297     end
0298 end
0299 if verbosityLevel
0300     if issparse(blockVectorX)
0301         fprintf(['The sparse initial guess with %i colunms '...
0302         'and %i raws is detected  \n'],n,blockSize);
0303     else
0304         fprintf(['The full initial guess with %i colunms '...
0305             'and %i raws is detected  \n'],n,blockSize);
0306     end
0307     if ischar(operatorA)
0308         fprintf('The main operator is detected as an M-function %s \n',...
0309             operatorA);
0310     elseif issparse(operatorA)
0311         fprintf('The main operator is detected as a sparse matrix \n');
0312     else
0313         fprintf('The main operator is detected as a full matrix \n');
0314     end
0315     if isempty(operatorB)
0316         fprintf('Solving standard eigenvalue problem, not generalized \n');
0317     elseif ischar(operatorB)
0318         fprintf(['The second operator of the generalized eigenproblem \n'...
0319         'is detected as an M-function %s \n'],operatorB);
0320     elseif issparse(operatorB)
0321         fprintf(strcat('The second operator of the generalized',... 
0322             'eigenproblem \n is detected as a sparse matrix \n'));
0323     else
0324         fprintf(strcat('The second operator of the generalized',... 
0325             'eigenproblem \n is detected as a full matrix \n'));        
0326     end
0327     if isempty(operatorT)
0328         fprintf('No preconditioner is detected \n');
0329     else
0330         fprintf('The preconditioner is detected as an M-function %s \n',...
0331             operatorT);
0332     end
0333     if isempty(blockVectorY)
0334         fprintf('No matrix of constraints is detected \n')
0335     elseif issparse(blockVectorY)
0336         fprintf('The sparse matrix of %i constraints is detected \n',sizeY);
0337     else
0338         fprintf('The full matrix of %i constraints is detected \n',sizeY);
0339     end
0340     if issparse(blockVectorY) ~= issparse(blockVectorX)
0341         warning('BLOPEX:lobpcg:SparsityInconsistent',...
0342             strcat('The sparsity formats of the initial guess and ',...
0343             'the constraints are inconsistent'));
0344     end
0345 end
0346 
0347 % Set defaults
0348 
0349 if isempty(residualTolerance)
0350     residualTolerance = sqrt(eps)*n;
0351 end
0352 if isempty(maxIterations)
0353     maxIterations = min(n,20);
0354 end
0355 if isempty(verbosityLevel)
0356     verbosityLevel = 0;
0357 end
0358 
0359 if verbosityLevel
0360     fprintf('Tolerance %e and maximum number of iterations %i \n',...
0361         residualTolerance,maxIterations)
0362 end
0363 
0364 %constraints preprocessing
0365 if isempty(blockVectorY)
0366     constraintStyle = 0;
0367 else
0368     %    constraintStyle = SYMMETRIC_CONSTRAINTS; % more accurate?
0369     constraintStyle = CONVENTIONAL_CONSTRAINTS;
0370 end
0371 
0372 if constraintStyle == CONVENTIONAL_CONSTRAINTS
0373     
0374     if isempty(operatorB)
0375         gramY = blockVectorY'*blockVectorY;
0376     else
0377         if ~ischar(operatorB)
0378             blockVectorBY = operatorB*blockVectorY;
0379         else
0380             blockVectorBY = feval(operatorB,blockVectorY);
0381         end
0382         gramY=blockVectorY'*blockVectorBY;
0383     end
0384     gramY=(gramY'+gramY)*0.5;
0385     if isempty(operatorB)
0386         blockVectorX = blockVectorX - ...
0387             blockVectorY*(gramY\(blockVectorY'*blockVectorX));
0388     else
0389         blockVectorX =blockVectorX - ...
0390             blockVectorY*(gramY\(blockVectorBY'*blockVectorX));
0391     end
0392     
0393 elseif constraintStyle == SYMMETRIC_CONSTRAINTS
0394     
0395     if ~isempty(operatorB)
0396         if ~ischar(operatorB)
0397             blockVectorY = operatorB*blockVectorY;
0398         else
0399             blockVectorY = feval(operatorB,blockVectorY);
0400         end
0401     end
0402     if isempty(operatorT)
0403         gramY = blockVectorY'*blockVectorY;
0404     else
0405         blockVectorTY = feval(operatorT,blockVectorY);
0406         gramY = blockVectorY'*blockVectorTY;
0407     end
0408     gramY=(gramY'+gramY)*0.5;
0409     if isempty(operatorT)
0410         blockVectorX = blockVectorX - ...
0411             blockVectorY*(gramY\(blockVectorY'*blockVectorX));
0412     else
0413         blockVectorX = blockVectorX - ...
0414             blockVectorTY*(gramY\(blockVectorY'*blockVectorX));
0415     end
0416     
0417 end
0418 
0419 %Making the initial vectors (operatorB-) orthonormal
0420 if isempty(operatorB)
0421     %[blockVectorX,gramXBX] = qr(blockVectorX,0);
0422     gramXBX=blockVectorX'*blockVectorX;
0423     if ~isreal(gramXBX)
0424         gramXBX=(gramXBX+gramXBX')*0.5;
0425     end
0426     [gramXBX,cholFlag]=chol(gramXBX);
0427     if  cholFlag ~= 0
0428         error('BLOPEX:lobpcg:ConstraintsTooTight',...
0429            'The initial approximation after constraints is not full rank');
0430     end
0431     blockVectorX = blockVectorX/gramXBX;
0432 else
0433     %[blockVectorX,blockVectorBX] = orth(operatorB,blockVectorX);
0434     if ~ischar(operatorB)
0435         blockVectorBX = operatorB*blockVectorX;
0436     else
0437         blockVectorBX = feval(operatorB,blockVectorX);
0438     end
0439     gramXBX=blockVectorX'*blockVectorBX;
0440     if ~isreal(gramXBX)
0441         gramXBX=(gramXBX+gramXBX')*0.5;
0442     end
0443     [gramXBX,cholFlag]=chol(gramXBX);
0444     if  cholFlag ~= 0
0445         error('BLOPEX:lobpcg:InitialNotFullRank',...
0446             sprintf('%s\n%s', ...
0447             'The initial approximation after constraints is not full rank',...
0448             'or/and operatorB is not positive definite'));
0449     end
0450     blockVectorX = blockVectorX/gramXBX;
0451     blockVectorBX = blockVectorBX/gramXBX;
0452 end
0453 
0454 % Checking if the problem is big enough for the algorithm,
0455 % i.e. n-sizeY > 5*blockSize
0456 % Theoretically, the algorithm should be able to run if
0457 % n-sizeY > 3*blockSize,
0458 % but the extreme cases might be unstable, so we use 5 instead of 3 here.
0459 if n-sizeY < 5*blockSize
0460     error('BLOPEX:lobpcg:MatrixTooSmall',...
0461         sprintf('%s\n%s', ...
0462     'The problem size is too small, relative to the block size.',... 
0463     'Try using eig() or eigs() instead.'));
0464 end
0465 
0466 % Preallocation
0467 residualNormsHistory=zeros(blockSize,maxIterations);
0468 lambdaHistory=zeros(blockSize,maxIterations+1);
0469 condestGhistory=zeros(1,maxIterations+1);
0470 
0471 blockVectorBR=zeros(n,blockSize);
0472 blockVectorAR=zeros(n,blockSize);
0473 blockVectorP=zeros(n,blockSize);
0474 blockVectorAP=zeros(n,blockSize);
0475 blockVectorBP=zeros(n,blockSize);
0476 
0477 %Initial settings for the loop
0478 %keyboard;
0479 if ischar(operatorA)
0480     blockVectorAX = operatorA*blockVectorX;
0481 else
0482     blockVectorAX = feval(operatorA,blockVectorX);
0483 end
0484 
0485 gramXAX = full(blockVectorX'*blockVectorAX);
0486 gramXAX = (gramXAX + gramXAX')*0.5;
0487 % eig(...,'chol') uses only the diagonal and upper triangle -
0488 % not true in MATLAB
0489 % Octave v3.2.3, eig() does not support inputting 'chol'
0490 [coordX,gramXAX]=eig(gramXAX,eye(blockSize));
0491 
0492 lambda=diag(gramXAX); %eig returms eigenvalues on the diagonal
0493 
0494 if issparse(blockVectorX)
0495     coordX=sparse(coordX);
0496 end
0497 
0498 blockVectorX  =  blockVectorX*coordX;
0499 blockVectorAX = blockVectorAX*coordX;
0500 if ~isempty(operatorB)
0501     blockVectorBX = blockVectorBX*coordX;
0502 end
0503 clear coordX
0504 
0505 condestGhistory(1)=-log10(eps)/2;  %if too small cause unnecessary restarts
0506 
0507 lambdaHistory(1:blockSize,1)=lambda;
0508 
0509 activeMask = true(blockSize,1);
0510 % currentBlockSize = blockSize; %iterate all
0511 %
0512 % restart=1;%steepest descent
0513 
0514 %The main part of the method is the loop of the CG method: begin
0515 for iterationNumber=1:maxIterations
0516     
0517     %     %Computing the active residuals
0518     %     if isempty(operatorB)
0519     %         if currentBlockSize > 1
0520     %             blockVectorR(:,activeMask)=blockVectorAX(:,activeMask) - ...
0521     %                 blockVectorX(:,activeMask)*spdiags(lambda(activeMask),0,currentBlockSize,currentBlockSize);
0522     %         else
0523     %             blockVectorR(:,activeMask)=blockVectorAX(:,activeMask) - ...
0524     %                 blockVectorX(:,activeMask)*lambda(activeMask);
0525     %                   %to make blockVectorR full when lambda is just a scalar
0526     %         end
0527     %     else
0528     %         if currentBlockSize > 1
0529     %             blockVectorR(:,activeMask)=blockVectorAX(:,activeMask) - ...
0530     %                 blockVectorBX(:,activeMask)*spdiags(lambda(activeMask),0,currentBlockSize,currentBlockSize);
0531     %         else
0532     %             blockVectorR(:,activeMask)=blockVectorAX(:,activeMask) - ...
0533     %                 blockVectorBX(:,activeMask)*lambda(activeMask);
0534     %                       %to make blockVectorR full when lambda is just a scalar
0535     %         end
0536     %     end
0537     
0538     %Computing all residuals
0539     if isempty(operatorB)
0540         if blockSize > 1
0541             blockVectorR = blockVectorAX - ...
0542                 blockVectorX*spdiags(lambda,0,blockSize,blockSize);
0543         else
0544             blockVectorR = blockVectorAX - blockVectorX*lambda;
0545             %to make blockVectorR full when lambda is just a scalar
0546         end
0547     else
0548         if blockSize > 1
0549             blockVectorR=blockVectorAX - ...
0550                 blockVectorBX*spdiags(lambda,0,blockSize,blockSize);
0551         else
0552             blockVectorR = blockVectorAX - blockVectorBX*lambda;
0553             %to make blockVectorR full when lambda is just a scalar
0554         end
0555     end
0556     
0557     %Satisfying the constraints for the active residulas
0558     if constraintStyle == SYMMETRIC_CONSTRAINTS
0559         if isempty(operatorT)
0560             blockVectorR(:,activeMask) = blockVectorR(:,activeMask) - ...
0561                 blockVectorY*(gramY\(blockVectorY'*...
0562                 blockVectorR(:,activeMask)));
0563         else
0564             blockVectorR(:,activeMask) = blockVectorR(:,activeMask) - ...
0565                 blockVectorY*(gramY\(blockVectorTY'*...
0566                 blockVectorR(:,activeMask)));
0567         end
0568     end
0569     
0570     residualNorms=full(sqrt(sum(conj(blockVectorR).*blockVectorR)'));
0571     residualNormsHistory(1:blockSize,iterationNumber)=residualNorms;
0572     
0573     %index antifreeze
0574     activeMask = full(residualNorms > residualTolerance) & activeMask;
0575     %activeMask = full(residualNorms > residualTolerance);
0576     %above allows vectors back into active, which causes problems with frosen Ps
0577     %activeMask = full(residualNorms > 0);      %iterate all, ignore freeze
0578     
0579     currentBlockSize = sum(activeMask);
0580     if  currentBlockSize == 0
0581         failureFlag=0; %all eigenpairs converged
0582         break
0583     end
0584     
0585     %Applying the preconditioner operatorT to the active residulas
0586     if ~isempty(operatorT)
0587         blockVectorR(:,activeMask) = ...
0588             feval(operatorT,blockVectorR(:,activeMask));
0589     end
0590     
0591     if constraintStyle == CONVENTIONAL_CONSTRAINTS
0592         if isempty(operatorB)
0593             blockVectorR(:,activeMask) = blockVectorR(:,activeMask) - ...
0594                 blockVectorY*(gramY\(blockVectorY'*...
0595                 blockVectorR(:,activeMask)));
0596         else
0597             blockVectorR(:,activeMask) = blockVectorR(:,activeMask) - ...
0598                 blockVectorY*(gramY\(blockVectorBY'*...
0599                 blockVectorR(:,activeMask)));
0600         end
0601     end
0602     
0603     %Making active (preconditioned) residuals orthogonal to blockVectorX
0604     if isempty(operatorB)
0605         blockVectorR(:,activeMask) = blockVectorR(:,activeMask) - ...
0606             blockVectorX*(blockVectorX'*blockVectorR(:,activeMask));
0607     else
0608         blockVectorR(:,activeMask) = blockVectorR(:,activeMask) - ...
0609             blockVectorX*(blockVectorBX'*blockVectorR(:,activeMask));
0610     end
0611     
0612     %Making active residuals orthonormal
0613     if isempty(operatorB)
0614         %[blockVectorR(:,activeMask),gramRBR]=...
0615         %qr(blockVectorR(:,activeMask),0); %to increase stability
0616         gramRBR=blockVectorR(:,activeMask)'*blockVectorR(:,activeMask);
0617         if ~isreal(gramRBR)
0618             gramRBR=(gramRBR+gramRBR')*0.5; 
0619         end
0620         [gramRBR,cholFlag]=chol(gramRBR);
0621         if  cholFlag == 0
0622             blockVectorR(:,activeMask) = blockVectorR(:,activeMask)/gramRBR;
0623         else
0624             warning('BLOPEX:lobpcg:ResidualNotFullRank',...
0625                 'The residual is not full rank.');
0626             break
0627         end
0628     else
0629         if ~ischar(operatorB)
0630             blockVectorBR(:,activeMask) = ...
0631                 operatorB*blockVectorR(:,activeMask);
0632         else
0633             blockVectorBR(:,activeMask) = ...
0634                 feval(operatorB,blockVectorR(:,activeMask));
0635         end
0636         gramRBR=blockVectorR(:,activeMask)'*blockVectorBR(:,activeMask);
0637         if ~isreal(gramRBR)
0638             gramRBR=(gramRBR+gramRBR')*0.5; 
0639         end
0640         [gramRBR,cholFlag]=chol(gramRBR);
0641         if  cholFlag == 0
0642             blockVectorR(:,activeMask) = ...
0643                 blockVectorR(:,activeMask)/gramRBR;
0644             blockVectorBR(:,activeMask) = ...
0645                 blockVectorBR(:,activeMask)/gramRBR;
0646         else
0647             warning('BLOPEX:lobpcg:ResidualNotFullRankOrElse',...
0648             strcat('The residual is not full rank or/and operatorB ',...
0649             'is not positive definite.'));
0650             break
0651         end
0652         
0653     end
0654     clear gramRBR;
0655     
0656     if ischar(operatorA)
0657         blockVectorAR(:,activeMask) = operatorA*blockVectorR(:,activeMask);
0658     else
0659         blockVectorAR(:,activeMask) = ...
0660             feval(operatorA,blockVectorR(:,activeMask));
0661     end
0662     
0663     if iterationNumber > 1
0664         
0665         %Making active conjugate directions orthonormal
0666         if isempty(operatorB)
0667             %[blockVectorP(:,activeMask),gramPBP] = qr(blockVectorP(:,activeMask),0);
0668             gramPBP=blockVectorP(:,activeMask)'*blockVectorP(:,activeMask);
0669             if ~isreal(gramPBP)
0670                 gramPBP=(gramPBP+gramPBP')*0.5; 
0671             end
0672             [gramPBP,cholFlag]=chol(gramPBP);
0673             if  cholFlag == 0
0674                 blockVectorP(:,activeMask) = ...
0675                     blockVectorP(:,activeMask)/gramPBP;
0676                 blockVectorAP(:,activeMask) = ...
0677                     blockVectorAP(:,activeMask)/gramPBP;
0678             else
0679                 warning('BLOPEX:lobpcg:DirectionNotFullRank',...
0680                     'The direction matrix is not full rank.');
0681                 break
0682             end
0683         else
0684             gramPBP=blockVectorP(:,activeMask)'*blockVectorBP(:,activeMask);
0685             if ~isreal(gramPBP)
0686                 gramPBP=(gramPBP+gramPBP')*0.5; 
0687             end
0688             [gramPBP,cholFlag]=chol(gramPBP);
0689             if  cholFlag == 0
0690                 blockVectorP(:,activeMask) = ...
0691                     blockVectorP(:,activeMask)/gramPBP;
0692                 blockVectorAP(:,activeMask) = ...
0693                     blockVectorAP(:,activeMask)/gramPBP;
0694                 blockVectorBP(:,activeMask) = ...
0695                     blockVectorBP(:,activeMask)/gramPBP;
0696             else
0697                 warning('BLOPEX:lobpcg:DirectionNotFullRank',...
0698                strcat('The direction matrix is not full rank ',...
0699             'or/and operatorB is not positive definite.'));
0700                 break
0701             end
0702         end
0703         clear gramPBP
0704     end
0705     
0706     condestGmean = mean(condestGhistory(max(1,iterationNumber-10-...
0707         round(log(currentBlockSize))):iterationNumber));
0708     
0709     %  restart=1;
0710     
0711     % The Raileight-Ritz method for [blockVectorX blockVectorR blockVectorP]
0712     
0713     if  residualNorms > eps^0.6
0714         explicitGramFlag = 0;
0715     else
0716         explicitGramFlag = 1;  %suggested by Garrett Moran, private
0717     end
0718     
0719     activeRSize=size(blockVectorR(:,activeMask),2);
0720     if iterationNumber == 1
0721         activePSize=0;
0722         restart=1;
0723     else
0724         activePSize=size(blockVectorP(:,activeMask),2);
0725         restart=0;
0726     end
0727     
0728     gramXAR=full(blockVectorAX'*blockVectorR(:,activeMask));
0729     gramRAR=full(blockVectorAR(:,activeMask)'*blockVectorR(:,activeMask));
0730     gramRAR=(gramRAR'+gramRAR)*0.5;
0731     
0732     if explicitGramFlag
0733         gramXAX=full(blockVectorAX'*blockVectorX);
0734         gramXAX=(gramXAX'+gramXAX)*0.5;
0735         if isempty(operatorB)
0736             gramXBX=full(blockVectorX'*blockVectorX);
0737             gramRBR=full(blockVectorR(:,activeMask)'*...
0738                 blockVectorR(:,activeMask));
0739             gramXBR=full(blockVectorX'*blockVectorR(:,activeMask));
0740         else
0741             gramXBX=full(blockVectorBX'*blockVectorX);
0742             gramRBR=full(blockVectorBR(:,activeMask)'*...
0743                 blockVectorR(:,activeMask));
0744             gramXBR=full(blockVectorBX'*blockVectorR(:,activeMask));
0745         end
0746         gramXBX=(gramXBX'+gramXBX)*0.5;
0747         gramRBR=(gramRBR'+gramRBR)*0.5;
0748         
0749     end
0750     
0751     for cond_try=1:2,           %cond_try == 2 when restart
0752         
0753         if ~restart
0754             gramXAP=full(blockVectorAX'*blockVectorP(:,activeMask));
0755             gramRAP=full(blockVectorAR(:,activeMask)'*...
0756                 blockVectorP(:,activeMask));
0757             gramPAP=full(blockVectorAP(:,activeMask)'*...
0758                 blockVectorP(:,activeMask));
0759             gramPAP=(gramPAP'+gramPAP)*0.5;
0760             
0761             if explicitGramFlag
0762                 gramA = [ gramXAX     gramXAR     gramXAP
0763                     gramXAR'    gramRAR     gramRAP
0764                     gramXAP'     gramRAP'    gramPAP ];
0765             else
0766                 gramA = [ diag(lambda)  gramXAR  gramXAP
0767                     gramXAR'      gramRAR  gramRAP
0768                     gramXAP'      gramRAP'  gramPAP ];
0769             end
0770             
0771             clear gramXAP  gramRAP gramPAP
0772             
0773             if isempty(operatorB)
0774                 gramXBP=full(blockVectorX'*blockVectorP(:,activeMask));
0775                 gramRBP=full(blockVectorR(:,activeMask)'*...
0776                     blockVectorP(:,activeMask));
0777             else
0778                 gramXBP=full(blockVectorBX'*blockVectorP(:,activeMask));
0779                 gramRBP=full(blockVectorBR(:,activeMask)'*...
0780                     blockVectorP(:,activeMask));
0781                 %or blockVectorR(:,activeMask)'*blockVectorBP(:,activeMask);
0782             end
0783             
0784             if explicitGramFlag
0785                 if isempty(operatorB)
0786                     gramPBP=full(blockVectorP(:,activeMask)'*...
0787                         blockVectorP(:,activeMask));
0788                 else
0789                     gramPBP=full(blockVectorBP(:,activeMask)'*...
0790                         blockVectorP(:,activeMask));
0791                 end
0792                 gramPBP=(gramPBP'+gramPBP)*0.5;
0793                 gramB = [ gramXBX  gramXBR  gramXBP
0794                     gramXBR' gramRBR  gramRBP
0795                     gramXBP' gramRBP' gramPBP ];
0796                 clear   gramPBP
0797             else
0798                 gramB=[eye(blockSize) zeros(blockSize,activeRSize) gramXBP
0799                     zeros(blockSize,activeRSize)' eye(activeRSize) gramRBP
0800                     gramXBP' gramRBP' eye(activePSize) ];
0801             end
0802             
0803             clear gramXBP  gramRBP;
0804             
0805         else
0806             
0807             if explicitGramFlag
0808                 gramA = [ gramXAX   gramXAR
0809                     gramXAR'    gramRAR  ];
0810                 gramB = [ gramXBX  gramXBR
0811                     gramXBR' eye(activeRSize)  ];
0812                 clear gramXAX gramXBX gramXBR
0813             else
0814                 gramA = [ diag(lambda)  gramXAR
0815                     gramXAR'        gramRAR  ];
0816                 gramB = eye(blockSize+activeRSize);
0817             end
0818             
0819             clear gramXAR gramRAR;
0820             
0821         end
0822         
0823         condestG = log10(cond(gramB))+1;
0824         if (condestG/condestGmean > 2 && condestG > 2 )|| condestG > 8
0825             %black magic - need to guess the restart
0826             if verbosityLevel
0827                 fprintf('Restart on step %i as condestG %5.4e \n',...
0828                     iterationNumber,condestG);
0829             end
0830             if cond_try == 1 && ~restart
0831                 restart=1; %steepest descent restart for stability
0832             else
0833                 warning('BLOPEX:lobpcg:IllConditioning',...
0834                     'Gramm matrix ill-conditioned: results unpredictable');
0835             end
0836         else
0837             break
0838         end
0839         
0840     end
0841     
0842     [gramA,gramB]=eig(gramA,gramB);
0843     lambda=diag(gramB(1:blockSize,1:blockSize)); 
0844     coordX=gramA(:,1:blockSize);
0845     
0846     clear gramA gramB
0847     
0848     if issparse(blockVectorX)
0849         coordX=sparse(coordX);
0850     end
0851     
0852     if ~restart
0853         blockVectorP =  blockVectorR(:,activeMask)*...
0854             coordX(blockSize+1:blockSize+activeRSize,:) + ...
0855             blockVectorP(:,activeMask)*...
0856             coordX(blockSize+activeRSize+1:blockSize + ...
0857             activeRSize+activePSize,:);
0858         blockVectorAP = blockVectorAR(:,activeMask)*...
0859             coordX(blockSize+1:blockSize+activeRSize,:) + ...
0860             blockVectorAP(:,activeMask)*...
0861             coordX(blockSize+activeRSize+1:blockSize + ...
0862             activeRSize+activePSize,:);
0863         if ~isempty(operatorB)
0864             blockVectorBP = blockVectorBR(:,activeMask)*...
0865                 coordX(blockSize+1:blockSize+activeRSize,:) + ...
0866                 blockVectorBP(:,activeMask)*...
0867                 coordX(blockSize+activeRSize+1:blockSize+activeRSize+activePSize,:);
0868         end
0869     else %use block steepest descent
0870         blockVectorP =   blockVectorR(:,activeMask)*...
0871             coordX(blockSize+1:blockSize+activeRSize,:);
0872         blockVectorAP = blockVectorAR(:,activeMask)*...
0873             coordX(blockSize+1:blockSize+activeRSize,:);
0874         if ~isempty(operatorB)
0875             blockVectorBP = blockVectorBR(:,activeMask)*...
0876                 coordX(blockSize+1:blockSize+activeRSize,:);
0877         end
0878     end
0879     
0880     blockVectorX = blockVectorX*coordX(1:blockSize,:) + blockVectorP;
0881     blockVectorAX=blockVectorAX*coordX(1:blockSize,:) + blockVectorAP;
0882     if ~isempty(operatorB)
0883         blockVectorBX=blockVectorBX*coordX(1:blockSize,:) + blockVectorBP;
0884     end
0885     clear coordX
0886     %%end RR
0887     
0888     lambdaHistory(1:blockSize,iterationNumber+1)=lambda;
0889     condestGhistory(iterationNumber+1)=condestG;
0890     
0891     if verbosityLevel
0892         fprintf('Iteration %i current block size %i \n',...
0893             iterationNumber,currentBlockSize);
0894         fprintf('Eigenvalues lambda %17.16e \n',lambda);
0895         fprintf('Residual Norms %e \n',residualNorms');
0896     end
0897 end
0898 %The main step of the method was the CG cycle: end
0899 
0900 %Postprocessing
0901 
0902 %Making sure blockVectorX's "exactly" satisfy the blockVectorY constrains??
0903 
0904 %Making sure blockVectorX's are "exactly" othonormalized by final "exact" RR
0905 if isempty(operatorB)
0906     gramXBX=full(blockVectorX'*blockVectorX);
0907 else
0908     if ~ischar(operatorB)
0909         blockVectorBX = operatorB*blockVectorX;
0910     else
0911         blockVectorBX = feval(operatorB,blockVectorX);
0912     end
0913     gramXBX=full(blockVectorX'*blockVectorBX);
0914 end
0915 gramXBX=(gramXBX'+gramXBX)*0.5;
0916 
0917 if ischar(operatorA)
0918     blockVectorAX = operatorA*blockVectorX;
0919 else
0920     blockVectorAX = feval(operatorA,blockVectorX);
0921 end
0922 gramXAX = full(blockVectorX'*blockVectorAX);
0923 gramXAX = (gramXAX + gramXAX')*0.5;
0924 
0925 %Raileigh-Ritz for blockVectorX, which is already operatorB-orthonormal
0926 [coordX,gramXBX] = eig(gramXAX,gramXBX);
0927 lambda=diag(gramXBX);
0928 
0929 if issparse(blockVectorX)
0930     coordX=sparse(coordX);
0931 end
0932 
0933 blockVectorX  =   blockVectorX*coordX;
0934 blockVectorAX  =  blockVectorAX*coordX;
0935 if ~isempty(operatorB)
0936     blockVectorBX  =  blockVectorBX*coordX;
0937 end
0938 
0939 %Computing all residuals
0940 if isempty(operatorB)
0941     if blockSize > 1
0942         blockVectorR = blockVectorAX - ...
0943             blockVectorX*spdiags(lambda,0,blockSize,blockSize);
0944     else
0945         blockVectorR = blockVectorAX - blockVectorX*lambda;
0946         %to make blockVectorR full when lambda is just a scalar
0947     end
0948 else
0949     if blockSize > 1
0950         blockVectorR=blockVectorAX - ...
0951             blockVectorBX*spdiags(lambda,0,blockSize,blockSize);
0952     else
0953         blockVectorR = blockVectorAX - blockVectorBX*lambda;
0954         %to make blockVectorR full when lambda is just a scalar
0955     end
0956 end
0957 
0958 residualNorms=full(sqrt(sum(conj(blockVectorR).*blockVectorR)'));
0959 residualNormsHistory(1:blockSize,iterationNumber)=residualNorms;
0960 
0961 if verbosityLevel
0962     fprintf('Final Eigenvalues lambda %17.16e \n',lambda);
0963     fprintf('Final Residual Norms %e \n',residualNorms');
0964 end
0965 
0966 if verbosityLevel == 2
0967     whos
0968     if matlabFlag==1;
0969         figure(491)
0970         semilogy((residualNormsHistory(1:blockSize,1:iterationNumber-1))');
0971         title('Residuals for Different Eigenpairs','fontsize',16);
0972         ylabel('Eucledian norm of residuals','fontsize',16);
0973         xlabel('Iteration number','fontsize',16);
0974         axis tight;
0975         %axis([0 maxIterations+1 1e-15 1e3])
0976         set(gca,'FontSize',14);
0977         figure(492);
0978         semilogy((lambdaHistory(1:blockSize,1:iterationNumber)-...
0979             repmat(lambda,1,iterationNumber))');
0980         title('Eigenvalue errors for Different Eigenpairs','fontsize',16);
0981         ylabel('Estimated eigenvalue errors','fontsize',16);
0982         xlabel('Iteration number','fontsize',16);
0983         axis tight;
0984         %axis([0 maxIterations+1 1e-15 1e3])
0985         set(gca,'FontSize',14);
0986         drawnow;
0987     end
0988 end
0989 
0990 varargout(1)={failureFlag};
0991 varargout(2)={lambdaHistory(1:blockSize,1:iterationNumber)};
0992 varargout(3)={residualNormsHistory(1:blockSize,1:iterationNumber-1)};
0993 return

Generated on Wed 08-Feb-2012 18:20:24 by m2html © 2005