static char help[] = "Computes the smallest nonzero eigenvalue of the Laplacian of a graph.\n\n"
"This example illustrates EPSAttachDeflationSpace(). The example graph corresponds to a "
"2-D regular mesh. The command line options are:\n\n"
" -n <n>, where <n> = number of grid subdivisions in x dimension.\n\n"
" -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";
#include "slepceps.h"
#undef __FUNCT__
#define __FUNCT__ "main"
int main( int argc, char **argv )
{
Mat A; /* operator matrix */
Vec x;
EPS eps; /* eigenproblem solver context */
ST st;
EPSType type;
PetscReal error, tol, re, im;
PetscScalar kr, ki;
int N, n=10, m, nev, ierr, maxit, i, j, I, J, its, nconv, Istart, Iend;
PetscScalar v, w, one = 1.0;
PetscTruth flag;
SlepcInitialize(&argc,&argv,(char*)0,help);
ierr = PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);CHKERRQ(ierr);
ierr = PetscOptionsGetInt(PETSC_NULL,"-m",&m,&flag);CHKERRQ(ierr);
if( flag==PETSC_FALSE ) m=n;
N = n*m;
ierr = PetscPrintf(PETSC_COMM_WORLD,"\nFiedler vector of a 2-D regular mesh, N=%d (%dx%d grid)\n\n",N,n,m);CHKERRQ(ierr);
/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
Compute the operator matrix that defines the eigensystem, Ax=kx
In this example, A = L(G), where L is the Laplacian of graph G, i.e.
Lii = degree of node i, Lij = -1 if edge (i,j) exists in G
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
ierr = MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&A);CHKERRQ(ierr);
ierr = MatSetFromOptions(A);CHKERRQ(ierr);
ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr);
for( I=Istart; I<Iend; I++ ) {
v = -1.0; i = I/n; j = I-i*n;
w = 0.0;
if(i>0) { J=I-n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr); w++; }
if(i<m-1) { J=I+n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr); w++; }
if(j>0) { J=I-1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr); w++; }
if(j<n-1) { J=I+1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr); w++; }
MatSetValues(A,1,&I,1,&I,&w,INSERT_VALUES);CHKERRQ(ierr);
}
ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);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 standard eigenvalue problem
*/
ierr = EPSSetOperators(eps,A,PETSC_NULL);CHKERRQ(ierr);
/*
Specify some options: use shift-and-invert to compute eigenpairs
close to the origin
*/
ierr = EPSGetST(eps,&st);CHKERRQ(ierr);
ierr = STSetType(st,STSINV);CHKERRQ(ierr);
ierr = STSetShift(st,0);CHKERRQ(ierr);
/*
Set solver parameters at runtime
*/
ierr = EPSSetFromOptions(eps);CHKERRQ(ierr);
/*
Attach deflation space: in this case, the matrix has a constant
nullspace, [1 1 ... 1]^T is the eigenvector of the zero eigenvalue
*/
ierr = MatGetVecs(A,&x,PETSC_NULL);CHKERRQ(ierr);
ierr = VecSet(&one,x);CHKERRQ(ierr);
ierr = EPSAttachDeflationSpace(eps,1,&x,PETSC_FALSE);CHKERRQ(ierr);
ierr = VecDestroy(x);
/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
Solve the eigensystem
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
ierr = EPSSolve(eps);CHKERRQ(ierr);
ierr = EPSGetIterationNumber(eps, &its);CHKERRQ(ierr);
ierr = PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);CHKERRQ(ierr);
/*
Optional: Get some information from the solver and display it
*/
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);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 approximate 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-kx||/||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);
#ifdef PETSC_USE_COMPLEX
re = PetscRealPart(kr);
im = PetscImaginaryPart(kr);
#else
re = kr;
im = ki;
#endif
if (im!=0.0) {
ierr = PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12f\n",re,im,error);CHKERRQ(ierr);
} else {
ierr = PetscPrintf(PETSC_COMM_WORLD," %12f %12f\n",re,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 = SlepcFinalize();CHKERRQ(ierr);
return 0;
}