Actual source code: ex2.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:       SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:       Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain

  6:       This file is part of SLEPc. See the README file for conditions of use
  7:       and additional information.
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 11: static char help[] = "Standard symmetric eigenproblem corresponding to the Laplacian operator in 2 dimensions.\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 14:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";

 16:  #include slepceps.h

 20: int main( int argc, char **argv )
 21: {
 22:   Mat                  A;                  /* operator matrix */
 23:   EPS                  eps;                  /* eigenproblem solver context */
 24:   EPSType              type;
 25:   PetscReal            error, tol, re, im;
 26:   PetscScalar          kr, ki;
 28:   PetscInt             N, n=10, m, Istart, Iend, II, J;
 29:   int                  nev, maxit, i, j, its, nconv;
 30:   PetscScalar          v;
 31:   PetscTruth           flag;

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

 35:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 36:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,&flag);
 37:   if( flag==PETSC_FALSE ) m=n;
 38:   N = n*m;
 39:   PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem, N=%d (%dx%d grid)\n\n",N,n,m);

 41:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 42:      Compute the operator matrix that defines the eigensystem, Ax=kx
 43:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 45:   MatCreate(PETSC_COMM_WORLD,&A);
 46:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 47:   MatSetFromOptions(A);
 48: 
 49:   MatGetOwnershipRange(A,&Istart,&Iend);
 50:   for( II=Istart; II<Iend; II++ ) {
 51:     v = -1.0; i = II/n; j = II-i*n;
 52:     if(i>0) { J=II-n; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES); }
 53:     if(i<m-1) { J=II+n; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES); }
 54:     if(j>0) { J=II-1; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES); }
 55:     if(j<n-1) { J=II+1; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES); }
 56:     v=4.0; MatSetValues(A,1,&II,1,&II,&v,INSERT_VALUES);
 57:   }

 59:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 60:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 62:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 63:                 Create the eigensolver and set various options
 64:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 66:   /* 
 67:      Create eigensolver context
 68:   */
 69:   EPSCreate(PETSC_COMM_WORLD,&eps);

 71:   /* 
 72:      Set operators. In this case, it is a standard eigenvalue problem
 73:   */
 74:   EPSSetOperators(eps,A,PETSC_NULL);
 75:   EPSSetProblemType(eps,EPS_HEP);

 77:   /*
 78:      Set solver parameters at runtime
 79:   */
 80:   EPSSetFromOptions(eps);

 82:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 83:                       Solve the eigensystem
 84:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 86:   EPSSolve(eps);
 87:   EPSGetIterationNumber(eps, &its);
 88:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);

 90:   /*
 91:      Optional: Get some information from the solver and display it
 92:   */
 93:   EPSGetType(eps,&type);
 94:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
 95:   EPSGetDimensions(eps,&nev,PETSC_NULL);
 96:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
 97:   EPSGetTolerances(eps,&tol,&maxit);
 98:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);

100:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
101:                     Display solution and clean up
102:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

104:   /* 
105:      Get number of converged approximate eigenpairs
106:   */
107:   EPSGetConverged(eps,&nconv);
108:   PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);
109: 

111:   if (nconv>0) {
112:     /*
113:        Display eigenvalues and relative errors
114:     */
115:     PetscPrintf(PETSC_COMM_WORLD,
116:          "           k          ||Ax-kx||/||kx||\n"
117:          "   ----------------- ------------------\n" );

119:     for( i=0; i<nconv; i++ ) {
120:       /* 
121:         Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
122:         ki (imaginary part)
123:       */
124:       EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
125:       /*
126:          Compute the relative error associated to each eigenpair
127:       */
128:       EPSComputeRelativeError(eps,i,&error);

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