Actual source code: ex1.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 1 dimension.\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = number of grid subdivisions = matrix dimension.\n\n";

 15:  #include slepceps.h

 19: int main( int argc, char **argv )
 20: {
 21:   Mat            A;           /* operator matrix */
 22:   EPS            eps;         /* eigenproblem solver context */
 23:   EPSType        type;
 24:   PetscReal      error, tol,re, im;
 25:   PetscScalar    kr, ki;
 26:   Vec            xr, xi;
 28:   PetscInt       n=30, i, Istart, Iend, col[3];
 29:   int            nev, maxit,its, nconv;
 30:   PetscTruth     FirstBlock=PETSC_FALSE, LastBlock=PETSC_FALSE;
 31:   PetscScalar    value[3];

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

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

 38:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 39:      Compute the operator matrix that defines the eigensystem, Ax=kx
 40:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 42:   MatCreate(PETSC_COMM_WORLD,&A);
 43:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 44:   MatSetFromOptions(A);
 45: 
 46:   MatGetOwnershipRange(A,&Istart,&Iend);
 47:   if (Istart==0) FirstBlock=PETSC_TRUE;
 48:   if (Iend==n) LastBlock=PETSC_TRUE;
 49:   value[0]=-1.0; value[1]=2.0; value[2]=-1.0;
 50:   for( i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++ ) {
 51:     col[0]=i-1; col[1]=i; col[2]=i+1;
 52:     MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 53:   }
 54:   if (LastBlock) {
 55:     i=n-1; col[0]=n-2; col[1]=n-1;
 56:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 57:   }
 58:   if (FirstBlock) {
 59:     i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0;
 60:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 61:   }

 63:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 64:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 66:   MatGetVecs(A,PETSC_NULL,&xr);
 67:   MatGetVecs(A,PETSC_NULL,&xi);

 69:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 70:                 Create the eigensolver and set various options
 71:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 72:   /* 
 73:      Create eigensolver context
 74:   */
 75:   EPSCreate(PETSC_COMM_WORLD,&eps);

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

 83:   /*
 84:      Set solver parameters at runtime
 85:   */
 86:   EPSSetFromOptions(eps);

 88:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 89:                       Solve the eigensystem
 90:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

105:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
106:                     Display solution and clean up
107:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108:   /* 
109:      Get number of converged approximate eigenpairs
110:   */
111:   EPSGetConverged(eps,&nconv);
112:   PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %d\n\n",nconv);

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

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

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