Actual source code: ex2.c

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

  6:    This file is part of SLEPc.
  7:       
  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

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

 27:  #include slepceps.h

 31: int main( int argc, char **argv )
 32: {
 33:   Mat                  A;                  /* operator matrix */
 34:   EPS                  eps;                  /* eigenproblem solver context */
 35:   const EPSType  type;
 36:   PetscReal            error, tol, re, im;
 37:   PetscScalar          kr, ki;
 39:   PetscInt             N, n=10, m, Istart, Iend, II, nev, maxit, i, j, its, nconv;
 40:   PetscTruth           flag;

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

 44:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 45:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,&flag);
 46:   if( flag==PETSC_FALSE ) m=n;
 47:   N = n*m;
 48:   PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem, N=%d (%dx%d grid)\n\n",N,n,m);

 50:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 51:      Compute the operator matrix that defines the eigensystem, Ax=kx
 52:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 54:   MatCreate(PETSC_COMM_WORLD,&A);
 55:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 56:   MatSetFromOptions(A);
 57: 
 58:   MatGetOwnershipRange(A,&Istart,&Iend);
 59:   for( II=Istart; II<Iend; II++ ) {
 60:     i = II/n; j = II-i*n;
 61:     if(i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
 62:     if(i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
 63:     if(j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
 64:     if(j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
 65:     MatSetValue(A,II,II,4.0,INSERT_VALUES);
 66:   }

 68:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 69:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 71:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 72:                 Create the eigensolver and set various options
 73:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 75:   /* 
 76:      Create eigensolver context
 77:   */
 78:   EPSCreate(PETSC_COMM_WORLD,&eps);

 80:   /* 
 81:      Set operators. In this case, it is a standard eigenvalue problem
 82:   */
 83:   EPSSetOperators(eps,A,PETSC_NULL);
 84:   EPSSetProblemType(eps,EPS_HEP);

 86:   /*
 87:      Set solver parameters at runtime
 88:   */
 89:   EPSSetFromOptions(eps);

 91:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 92:                       Solve the eigensystem
 93:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 95:   EPSSolve(eps);
 96:   EPSGetIterationNumber(eps, &its);
 97:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);

 99:   /*
100:      Optional: Get some information from the solver and display it
101:   */
102:   EPSGetType(eps,&type);
103:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
104:   EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);
105:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
106:   EPSGetTolerances(eps,&tol,&maxit);
107:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);

109:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
110:                     Display solution and clean up
111:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

113:   /* 
114:      Get number of converged approximate eigenpairs
115:   */
116:   EPSGetConverged(eps,&nconv);
117:   PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);
118: 

120:   if (nconv>0) {
121:     /*
122:        Display eigenvalues and relative errors
123:     */
124:     PetscPrintf(PETSC_COMM_WORLD,
125:          "           k          ||Ax-kx||/||kx||\n"
126:          "   ----------------- ------------------\n" );

128:     for( i=0; i<nconv; i++ ) {
129:       /* 
130:         Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
131:         ki (imaginary part)
132:       */
133:       EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
134:       /*
135:          Compute the relative error associated to each eigenpair
136:       */
137:       EPSComputeRelativeError(eps,i,&error);

139: #ifdef PETSC_USE_COMPLEX
140:       re = PetscRealPart(kr);
141:       im = PetscImaginaryPart(kr);
142: #else
143:       re = kr;
144:       im = ki;
145: #endif 
146:       if (im!=0.0) {
147:         PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12g\n",re,im,error);
148:       } else {
149:         PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12g\n",re,error);
150:       }
151:     }
152:     PetscPrintf(PETSC_COMM_WORLD,"\n" );
153:   }
154: 
155:   /* 
156:      Free work space
157:   */
158:   EPSDestroy(eps);
159:   MatDestroy(A);
160:   SlepcFinalize();
161:   return 0;
162: }