Actual source code: ex2.c

  2: static char help[] = "Standard symmetric eigenproblem corresponding to the Laplacian operator in 2 dimensions.\n\n"
  3:   "The command line options are:\n"
  4:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\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;
 19:   PetscInt             N, n=10, m, Istart, Iend, II, J;
 20:   int                  nev, maxit, i, j, its, nconv;
 21:   PetscScalar          v;
 22:   PetscTruth           flag;

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

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

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

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

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

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

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

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

 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);

 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:   EPSGetTolerances(eps,&tol,&maxit);
 89:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);

 91:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 92:                     Display solution and clean up
 93:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

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