Actual source code: test6.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[] = "Diagonal eigenproblem.\n\n"
 23:   "The command line options are:\n"
 24:   "  -n <n>, where <n> = number of grid subdivisions = matrix dimension.\n"
 25:   "  -seed <s>, where <s> = seed for random number generation.\n\n";

 27: #include <slepceps.h>

 31: int main(int argc,char **argv)
 32: {
 33:   Mat            A;           /* problem matrix */
 34:   EPS            eps;         /* eigenproblem solver context */
 35:   Vec            v0;          /* initial vector */
 36:   PetscRandom    rand;
 37:   PetscReal      tol=1000*PETSC_MACHINE_EPSILON;
 38:   PetscInt       n=30,i,Istart,Iend,seed=0x12345678;

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

 43:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 44:   PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Eigenproblem, n=%d\n\n",n);

 46:   MatCreate(PETSC_COMM_WORLD,&A);
 47:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 48:   MatSetFromOptions(A);
 49:   MatGetOwnershipRange(A,&Istart,&Iend);
 50:   for (i=Istart;i<Iend;i++) {
 51:     MatSetValue(A,i,i,i+1,INSERT_VALUES);
 52:   }
 53:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 54:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 56:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 57:                       Solve the eigensystem
 58:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 59:   EPSCreate(PETSC_COMM_WORLD,&eps);
 60:   EPSSetOperators(eps,A,PETSC_NULL);
 61:   EPSSetProblemType(eps,EPS_HEP);
 62:   EPSSetTolerances(eps,tol,PETSC_DECIDE);
 63:   EPSSetFromOptions(eps);
 64:   /* set random initial vector */
 65:   MatGetVecs(A,&v0,PETSC_NULL);
 66:   PetscRandomCreate(PETSC_COMM_WORLD,&rand);
 67:   PetscRandomSetFromOptions(rand);
 68:   PetscOptionsGetInt(PETSC_NULL,"-seed",&seed,PETSC_NULL);
 69:   PetscRandomSetSeed(rand,seed);
 70:   PetscRandomSeed(rand);
 71:   VecSetRandom(v0,rand);
 72:   EPSSetInitialSpace(eps,1,&v0);
 73:   /* call the solver */
 74:   EPSSolve(eps);

 76:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 77:                     Display solution and clean up
 78:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 79:   EPSPrintSolution(eps,PETSC_NULL);
 80:   EPSDestroy(&eps);
 81:   MatDestroy(&A);
 82:   VecDestroy(&v0);
 83:   PetscRandomDestroy(&rand);
 84:   SlepcFinalize();
 85:   return 0;
 86: }