Actual source code: ex1.c

  2: static char help[] = "Solves a standard symmetric eigenproblem corresponding to the Laplacian operator in 1 dimension.\n\n"
  3:   "The command line options are:\n\n"
  4:   "  -n <n>, where <n> = number of grid subdivisions = matrix dimension.\n\n";

 6:  #include slepceps.h

 10: int main( int argc, char **argv )
 11: {
 12:   Mat         A;               /* operator matrix */
 13:   EPS         eps;             /* eigenproblem solver context */
 14:   EPSType     type;
 15:   PetscReal   error, tol,re, im;
 16:   PetscScalar kr, ki;
 17:   int         n=30, nev, ierr, maxit, i, its, nconv,
 18:               col[3], Istart, Iend, FirstBlock=0, LastBlock=0;
 19:   PetscScalar value[3];

 21:   SlepcInitialize(&argc,&argv,(char*)0,help);

 23:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 24:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian Eigenproblem, n=%d\n\n",n);
 25: 

 27:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 28:      Compute the operator matrix that defines the eigensystem, Ax=kx
 29:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 31:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,&A);
 32:   MatSetFromOptions(A);
 33: 
 34:   MatGetOwnershipRange(A,&Istart,&Iend);
 35:   if (Istart==0) FirstBlock=PETSC_TRUE;
 36:   if (Iend==n) LastBlock=PETSC_TRUE;
 37:   value[0]=-1.0; value[1]=2.0; value[2]=-1.0;
 38:   for( i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++ ) {
 39:     col[0]=i-1; col[1]=i; col[2]=i+1;
 40:     MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 41:   }
 42:   if (LastBlock) {
 43:     i=n-1; col[0]=n-2; col[1]=n-1;
 44:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 45:   }
 46:   if (FirstBlock) {
 47:     i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0;
 48:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 49:   }

 51:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 52:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 54:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 55:                 Create the eigensolver and set various options
 56:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 58:   /* 
 59:      Create eigensolver context
 60:   */
 61:   EPSCreate(PETSC_COMM_WORLD,&eps);

 63:   /* 
 64:      Set operators. In this case, it is a standard eigenvalue problem
 65:   */
 66:   EPSSetOperators(eps,A,PETSC_NULL);

 68:   /*
 69:      Set solver parameters at runtime
 70:   */
 71:   EPSSetFromOptions(eps);

 73:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 74:                       Solve the eigensystem
 75:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 77:   EPSSolve(eps);
 78:   EPSGetIterationNumber(eps, &its);
 79:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
 80: 
 81:   /*
 82:      Optional: Get some information from the solver and display it
 83:   */
 84:   EPSGetType(eps,&type);
 85:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
 86:   EPSGetDimensions(eps,&nev,PETSC_NULL);
 87:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
 88: 
 89:   EPSGetTolerances(eps,&tol,&maxit);
 90:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
 91: 


 94:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 95:                     Display solution and clean up
 96:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 98:   /* 
 99:      Get number of converged approximate eigenpairs
100:   */
101:   EPSGetConverged(eps,&nconv);
102:   PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %d\n\n",nconv);
103: 

105:   if (nconv>0) {
106:     /*
107:        Display eigenvalues and relative errors
108:     */
109:     PetscPrintf(PETSC_COMM_WORLD,
110:          "           k          ||Ax-kx||/||kx||\n"
111:          "   ----------------- ------------------\n" );

113:     for( i=0; i<nconv; i++ ) {
114:       /* 
115:         Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
116:         ki (imaginary part)
117:       */
118:       EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
119:       /*
120:          Compute the relative error associated to each eigenpair
121:       */
122:       EPSComputeRelativeError(eps,i,&error);

124: #ifdef PETSC_USE_COMPLEX
125:       re = PetscRealPart(kr);
126:       im = PetscImaginaryPart(kr);
127: #else
128:       re = kr;
129:       im = ki;
130: #endif 
131:       if (im!=0.0) {
132:         PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12f\n",re,im,error);
133:       } else {
134:         PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12f\n",re,error);
135:       }
136:     }
137:     PetscPrintf(PETSC_COMM_WORLD,"\n" );
138:   }
139: 
140:   /* 
141:      Free work space
142:   */
143:   EPSDestroy(eps);
144:   MatDestroy(A);
145:   SlepcFinalize();
146:   return 0;
147: }