/*
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
SLEPc - Scalable Library for Eigenvalue Problem Computations
Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain
This file is part of SLEPc.
SLEPc is free software: you can redistribute it and/or modify it under the
terms of version 3 of the GNU Lesser General Public License as published by
the Free Software Foundation.
SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
more details.
You should have received a copy of the GNU Lesser General Public License
along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
*/
static char help[] = "Solves a generalized eigensystem Ax=kBx with matrices loaded from a file.\n"
"This example works for both real and complex numbers.\n\n"
"The command line options are:\n"
" -f1 <filename>, where <filename> = matrix (A) file in PETSc binary form.\n"
" -f2 <filename>, where <filename> = matrix (B) file in PETSc binary form.\n\n";
#include "slepceps.h"
#undef __FUNCT__
#define __FUNCT__ "main"
int main( int argc, char **argv )
{
Mat A,B; /* matrices */
EPS eps; /* eigenproblem solver context */
const EPSType type;
PetscReal error, tol, re, im;
PetscScalar kr, ki;
PetscErrorCode ierr;
PetscInt nev, maxit, i, its, lits, nconv;
char filename[256];
PetscViewer viewer;
PetscTruth flg;
SlepcInitialize(&argc,&argv,(char*)0,help);
/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
Load the matrices that define the eigensystem, Ax=kBx
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
ierr = PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized eigenproblem stored in file.\n\n");CHKERRQ(ierr);
ierr = PetscOptionsGetString(PETSC_NULL,"-f1",filename,256,&flg);CHKERRQ(ierr);
if (!flg) {
SETERRQ(1,"Must indicate a file name for matrix A with the -f1 option.");
}
#if defined(PETSC_USE_COMPLEX)
ierr = PetscPrintf(PETSC_COMM_WORLD," Reading COMPLEX matrices from binary files...\n");CHKERRQ(ierr);
#else
ierr = PetscPrintf(PETSC_COMM_WORLD," Reading REAL matrices from binary files...\n");CHKERRQ(ierr);
#endif
ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
ierr = MatLoad(viewer,MATAIJ,&A);CHKERRQ(ierr);
ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
ierr = PetscOptionsGetString(PETSC_NULL,"-f2",filename,256,&flg);CHKERRQ(ierr);
if (!flg) {
SETERRQ(1,"Must indicate a file name for matrix B with the -f2 option.");
}
ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
ierr = MatLoad(viewer,MATAIJ,&B);CHKERRQ(ierr);
ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
Create the eigensolver and set various options
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
/*
Create eigensolver context
*/
ierr = EPSCreate(PETSC_COMM_WORLD,&eps);CHKERRQ(ierr);
/*
Set operators. In this case, it is a generalized eigenvalue problem
*/
ierr = EPSSetOperators(eps,A,B);CHKERRQ(ierr);
/*
Set solver parameters at runtime
*/
ierr = EPSSetFromOptions(eps);CHKERRQ(ierr);
/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
Solve the eigensystem
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
ierr = EPSSolve(eps);CHKERRQ(ierr);
/*
Optional: Get some information from the solver and display it
*/
ierr = EPSGetIterationNumber(eps, &its);CHKERRQ(ierr);
ierr = PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);CHKERRQ(ierr);
ierr = EPSGetOperationCounters(eps,PETSC_NULL,PETSC_NULL,&lits);CHKERRQ(ierr);
ierr = PetscPrintf(PETSC_COMM_WORLD," Number of linear iterations of the method: %d\n",lits);CHKERRQ(ierr);
ierr = EPSGetType(eps,&type);CHKERRQ(ierr);
ierr = PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);CHKERRQ(ierr);
ierr = EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
ierr = PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);CHKERRQ(ierr);
ierr = EPSGetTolerances(eps,&tol,&maxit);CHKERRQ(ierr);
ierr = PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);CHKERRQ(ierr);
/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
Display solution and clean up
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
/*
Get number of converged eigenpairs
*/
ierr = EPSGetConverged(eps,&nconv);CHKERRQ(ierr);
ierr = PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);CHKERRQ(ierr);
if (nconv>0) {
/*
Display eigenvalues and relative errors
*/
ierr = PetscPrintf(PETSC_COMM_WORLD,
" k ||Ax-kBx||/||kx||\n"
" --------------------- ------------------\n" );CHKERRQ(ierr);
for( i=0; i<nconv; i++ ) {
/*
Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
ki (imaginary part)
*/
ierr = EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
/*
Compute the relative error associated to each eigenpair
*/
ierr = EPSComputeRelativeError(eps,i,&error);CHKERRQ(ierr);
#if defined(PETSC_USE_COMPLEX)
re = PetscRealPart(kr);
im = PetscImaginaryPart(kr);
#else
re = kr;
im = ki;
#endif
if( im != 0.0 ) {
ierr = PetscPrintf(PETSC_COMM_WORLD," % 6f %+6f i",re,im);CHKERRQ(ierr);
} else {
ierr = PetscPrintf(PETSC_COMM_WORLD," % 6f ",re); CHKERRQ(ierr);
}
ierr = PetscPrintf(PETSC_COMM_WORLD," % 12g\n",error);CHKERRQ(ierr);
}
ierr = PetscPrintf(PETSC_COMM_WORLD,"\n" );CHKERRQ(ierr);
}
/*
Free work space
*/
ierr = EPSDestroy(eps);CHKERRQ(ierr);
ierr = MatDestroy(A);CHKERRQ(ierr);
ierr = MatDestroy(B);CHKERRQ(ierr);
ierr = SlepcFinalize();CHKERRQ(ierr);
return 0;
}