Actual source code: ex2.c

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

 7:  #include slepceps.h

 11: int main( int argc, char **argv )
 12: {
 13:   Mat         A;               /* operator matrix */
 14:   EPS         eps;             /* eigenproblem solver context */
 15:   EPSType     type;
 16:   PetscReal   error, tol, re, im;
 17:   PetscScalar kr, ki;
 18:   int         N, n=10, m, nev, ierr, maxit, i, j, I, J, its, nconv, Istart, Iend;
 19:   PetscScalar v;
 20:   PetscTruth  flag;

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

 24:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 25:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,&flag);
 26:   if( flag==PETSC_FALSE ) m=n;
 27:   N = n*m;
 28:   PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem, N=%d (%dx%d grid)\n\n",N,n,m);

 30:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 31:      Compute the operator matrix that defines the eigensystem, Ax=kx
 32:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 34:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&A);
 35:   MatSetFromOptions(A);
 36: 
 37:   MatGetOwnershipRange(A,&Istart,&Iend);
 38:   for( I=Istart; I<Iend; I++ ) {
 39:     v = -1.0; i = I/n; j = I-i*n;
 40:     if(i>0) { J=I-n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES); }
 41:     if(i<m-1) { J=I+n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES); }
 42:     if(j>0) { J=I-1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES); }
 43:     if(j<n-1) { J=I+1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES); }
 44:     v=4.0; MatSetValues(A,1,&I,1,&I,&v,INSERT_VALUES);
 45:   }

 47:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 48:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 50:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 51:                 Create the eigensolver and set various options
 52:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 54:   /* 
 55:      Create eigensolver context
 56:   */
 57:   EPSCreate(PETSC_COMM_WORLD,&eps);

 59:   /* 
 60:      Set operators. In this case, it is a standard eigenvalue problem
 61:   */
 62:   EPSSetOperators(eps,A,PETSC_NULL);
 63:   EPSSetProblemType(eps,EPS_HEP);

 65:   /*
 66:      Set solver parameters at runtime
 67:   */
 68:   EPSSetFromOptions(eps);

 70:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 71:                       Solve the eigensystem
 72:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 74:   EPSSolve(eps);
 75:   EPSGetIterationNumber(eps, &its);
 76:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);

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

 88:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 89:                     Display solution and clean up
 90:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 92:   /* 
 93:      Get number of converged approximate eigenpairs
 94:   */
 95:   EPSGetConverged(eps,&nconv);
 96:   PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);
 97: 

 99:   if (nconv>0) {
100:     /*
101:        Display eigenvalues and relative errors
102:     */
103:     PetscPrintf(PETSC_COMM_WORLD,
104:          "           k          ||Ax-kx||/||kx||\n"
105:          "   ----------------- ------------------\n" );

107:     for( i=0; i<nconv; i++ ) {
108:       /* 
109:         Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
110:         ki (imaginary part)
111:       */
112:       EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
113:       /*
114:          Compute the relative error associated to each eigenpair
115:       */
116:       EPSComputeRelativeError(eps,i,&error);

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