Actual source code: ex3.c

  2: static char help[] = "Solves the same eigenproblem as in example ex2, but using a shell matrix. "
  3:   "The problem is a standard symmetric eigenproblem corresponding to the 2-D Laplacian operator.\n\n"
  4:   "The command line options are:\n\n"
  5:   "  -n <n>, where <n> = number of grid subdivisions in both x and y dimensions.\n\n";

 7:  #include slepceps.h
  8: #include "petscblaslapack.h"

 10: /* 
 11:    User-defined routines
 12: */
 13: extern int MatLaplacian2D_Mult( Mat A, Vec x, Vec y );

 17: int main( int argc, char **argv )
 18: {
 19:   Mat         A;               /* operator matrix */
 20:   EPS         eps;             /* eigenproblem solver context */
 21:   EPSType     type;
 22:   PetscReal   error, tol, re, im;
 23:   PetscScalar kr, ki;
 24:   int         size, N, n=10, nev, ierr, maxit, i, its, nconv;

 26:   SlepcInitialize(&argc,&argv,(char*)0,help);
 27:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 28:   if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");

 30:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 31:   N = n*n;
 32:   PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem (matrix-free version), N=%d (%dx%d grid)\n\n",N,n,n);

 34:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 35:      Compute the operator matrix that defines the eigensystem, Ax=kx
 36:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 38:   MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,&n,&A);
 39:   MatSetFromOptions(A);
 40:   MatShellSetOperation(A,MATOP_MULT,(void(*)())MatLaplacian2D_Mult);
 41:   MatShellSetOperation(A,MATOP_MULT_TRANSPOSE,(void(*)())MatLaplacian2D_Mult);

 43:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 44:                 Create the eigensolver and set various options
 45:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 47:   /* 
 48:      Create eigensolver context
 49:   */
 50:   EPSCreate(PETSC_COMM_WORLD,&eps);

 52:   /* 
 53:      Set operators. In this case, it is a standard eigenvalue problem
 54:   */
 55:   EPSSetOperators(eps,A,PETSC_NULL);

 57:   /*
 58:      Set solver parameters at runtime
 59:   */
 60:   EPSSetFromOptions(eps);

 62:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 63:                       Solve the eigensystem
 64:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 66:   EPSSolve(eps);
 67:   EPSGetIterationNumber(eps, &its);
 68:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);

 70:   /*
 71:      Optional: Get some information from the solver and display it
 72:   */
 73:   EPSGetType(eps,&type);
 74:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
 75:   EPSGetDimensions(eps,&nev,PETSC_NULL);
 76:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
 77:   EPSGetTolerances(eps,&tol,&maxit);
 78:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);

 80:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 81:                     Display solution and clean up
 82:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 84:   /* 
 85:      Get number of converged approximate eigenpairs
 86:   */
 87:   EPSGetConverged(eps,&nconv);
 88:   PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %d\n\n",nconv);
 89: 

 91:   if (nconv>0) {
 92:     /*
 93:        Display eigenvalues and relative errors
 94:     */
 95:     PetscPrintf(PETSC_COMM_WORLD,
 96:          "           k          ||Ax-kx||/||kx||\n"
 97:          "   ----------------- ------------------\n" );

 99:     for( i=0; i<nconv; i++ ) {
100:       /* 
101:         Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
102:         ki (imaginary part)
103:       */
104:       EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
105:       /*
106:          Compute the relative error associated to each eigenpair
107:       */
108:       EPSComputeRelativeError(eps,i,&error);

110: #ifdef PETSC_USE_COMPLEX
111:       re = PetscRealPart(kr);
112:       im = PetscImaginaryPart(kr);
113: #else
114:       re = kr;
115:       im = ki;
116: #endif 
117:       if (im!=0.0) {
118:         PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12f\n",re,im,error);
119:       } else {
120:         PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12f\n",re,error);
121:       }
122:     }
123:     PetscPrintf(PETSC_COMM_WORLD,"\n" );
124:   }
125: 
126:   /* 
127:      Free work space
128:   */
129:   EPSDestroy(eps);
130:   MatDestroy(A);
131:   SlepcFinalize();
132:   return 0;
133: }

135: /*
136:     Compute the matrix vector multiplication y<---T*x where T is a nx by nx
137:     tridiagonal matrix with DD on the diagonal, DL on the subdiagonal, and 
138:     DU on the superdiagonal.
139:  */
140: static void tv( int nx, PetscScalar *x, PetscScalar *y )
141: {
142:   PetscScalar dd, dl, du;
143:   int         j;

145:   dd  = 4.0;
146:   dl  = -1.0;
147:   du  = -1.0;

149:   y[0] =  dd*x[0] + du*x[1];
150:   for( j=1; j<nx-1; j++ )
151:     y[j] = dl*x[j-1] + dd*x[j] + du*x[j+1];
152:   y[nx-1] = dl*x[nx-2] + dd*x[nx-1];
153: }

157: /*
158:     Matrix-vector product subroutine for the 2D Laplacian.

160:     The matrix used is the 2 dimensional discrete Laplacian on unit square with
161:     zero Dirichlet boundary condition.
162:  
163:     Computes y <-- A*x, where A is the block tridiagonal matrix
164:  
165:                  | T -I          | 
166:                  |-I  T -I       |
167:              A = |   -I  T       |
168:                  |        ...  -I|
169:                  |           -I T|
170:  
171:     The subroutine TV is called to compute y<--T*x.
172:  */
173: int MatLaplacian2D_Mult( Mat A, Vec x, Vec y )
174: {
175:   void        *ctx;
176:   int         ierr, nx, lo, j, one=1;
177:   PetscScalar *px, *py, dmone=-1.0;
178: 
179:   MatShellGetContext( A, &ctx );
180:   nx = *(int *)ctx;
181:   VecGetArray( x, &px );
182:   VecGetArray( y, &py );

184:   tv( nx, &px[0], &py[0] );
185:   BLaxpy_( &nx, &dmone, &px[nx], &one, &py[0], &one );

187:   for( j=2; j<nx; j++ ) {
188:     lo = (j-1)*nx;
189:     tv( nx, &px[lo], &py[lo]);
190:     BLaxpy_( &nx, &dmone, &px[lo-nx], &one, &py[lo], &one );
191:     BLaxpy_( &nx, &dmone, &px[lo+nx], &one, &py[lo], &one );
192:   }

194:   lo = (nx-1)*nx;
195:   tv( nx, &px[lo], &py[lo]);
196:   BLaxpy_( &nx, &dmone, &px[lo-nx], &one, &py[lo], &one );

198:   VecRestoreArray( x, &px );
199:   VecRestoreArray( y, &py );

201:   return 0;
202: }