Actual source code: ex1.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2011, Universitat 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 1 dimension.\n\n"
 23:   "The command line options are:\n"
 24:   "  -n <n>, where <n> = number of grid subdivisions = matrix dimension.\n\n";

 26: #include <slepceps.h>

 30: int main(int argc,char **argv)
 31: {
 32:   Mat            A;           /* problem matrix */
 33:   EPS            eps;         /* eigenproblem solver context */
 34:   const EPSType  type;
 35:   PetscReal      error,tol,re,im;
 36:   PetscScalar    kr,ki,value[3];
 37:   Vec            xr,xi;
 38:   PetscInt       n=30,i,Istart,Iend,col[3],nev,maxit,its,nconv;
 39:   PetscBool      FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;

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

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

 47:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 48:      Compute the operator matrix that defines the eigensystem, Ax=kx
 49:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 51:   MatCreate(PETSC_COMM_WORLD,&A);
 52:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 53:   MatSetFromOptions(A);
 54: 
 55:   MatGetOwnershipRange(A,&Istart,&Iend);
 56:   if (Istart==0) FirstBlock=PETSC_TRUE;
 57:   if (Iend==n) LastBlock=PETSC_TRUE;
 58:   value[0]=-1.0; value[1]=2.0; value[2]=-1.0;
 59:   for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) {
 60:     col[0]=i-1; col[1]=i; col[2]=i+1;
 61:     MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 62:   }
 63:   if (LastBlock) {
 64:     i=n-1; col[0]=n-2; col[1]=n-1;
 65:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 66:   }
 67:   if (FirstBlock) {
 68:     i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0;
 69:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 70:   }

 72:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 73:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 75:   MatGetVecs(A,PETSC_NULL,&xr);
 76:   MatGetVecs(A,PETSC_NULL,&xi);

 78:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 79:                 Create the eigensolver and set various options
 80:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 81:   /* 
 82:      Create eigensolver context
 83:   */
 84:   EPSCreate(PETSC_COMM_WORLD,&eps);

 86:   /* 
 87:      Set operators. In this case, it is a standard eigenvalue problem
 88:   */
 89:   EPSSetOperators(eps,A,PETSC_NULL);
 90:   EPSSetProblemType(eps,EPS_HEP);

 92:   /*
 93:      Set solver parameters at runtime
 94:   */
 95:   EPSSetFromOptions(eps);

 97:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 98:                       Solve the eigensystem
 99:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

101:   EPSSolve(eps);
102:   /*
103:      Optional: Get some information from the solver and display it
104:   */
105:   EPSGetIterationNumber(eps,&its);
106:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",its);
107:   EPSGetType(eps,&type);
108:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
109:   EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);
110:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
111:   EPSGetTolerances(eps,&tol,&maxit);
112:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4G, maxit=%D\n",tol,maxit);

114:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
115:                     Display solution and clean up
116:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
117:   /* 
118:      Get number of converged approximate eigenpairs
119:   */
120:   EPSGetConverged(eps,&nconv);
121:   PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %D\n\n",nconv);

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

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

142: #if defined(PETSC_USE_COMPLEX)
143:       re = PetscRealPart(kr);
144:       im = PetscImaginaryPart(kr);
145: #else
146:       re = kr;
147:       im = ki;
148: #endif 
149:       if (im!=0.0) {
150:         PetscPrintf(PETSC_COMM_WORLD," %9F%+9F j %12G\n",re,im,error);
151:       } else {
152:         PetscPrintf(PETSC_COMM_WORLD,"   %12F       %12G\n",re,error);
153:       }
154:     }
155:     PetscPrintf(PETSC_COMM_WORLD,"\n");
156:   }
157: 
158:   /* 
159:      Free work space
160:   */
161:   EPSDestroy(&eps);
162:   MatDestroy(&A);
163:   VecDestroy(&xr);
164:   VecDestroy(&xi);
165:   SlepcFinalize();
166:   return 0;
167: }