Actual source code: test12.c

slepc-3.7.3 2016-09-29
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  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. Illustrates use of shell preconditioner.\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: PetscErrorCode PCApply_User(PC pc,Vec x,Vec y)
 32: {

 36:   VecCopy(x,y);
 37:   return(0);
 38: }

 42: int main(int argc,char **argv)
 43: {
 44:   Mat            A;           /* problem matrix */
 45:   EPS            eps;         /* eigenproblem solver context */
 46:   Vec            v0;          /* initial vector */
 47:   PetscRandom    rand;
 48:   PetscReal      tol=1000*PETSC_MACHINE_EPSILON;
 49:   PetscInt       n=30,i,Istart,Iend,seed=0x12345678;
 51:   ST             st;
 52:   KSP            ksp;
 53:   PC             pc;

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

 57:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 58:   PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Eigenproblem, n=%D\n\n",n);

 60:   MatCreate(PETSC_COMM_WORLD,&A);
 61:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 62:   MatSetFromOptions(A);
 63:   MatSetUp(A);
 64:   MatGetOwnershipRange(A,&Istart,&Iend);
 65:   for (i=Istart;i<Iend;i++) {
 66:     MatSetValue(A,i,i,i+1,INSERT_VALUES);
 67:   }
 68:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 69:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 71:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 72:                       Solve the eigensystem
 73:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 74:   EPSCreate(PETSC_COMM_WORLD,&eps);
 75:   EPSSetOperators(eps,A,NULL);
 76:   EPSSetProblemType(eps,EPS_HEP);
 77:   EPSSetTolerances(eps,tol,PETSC_DEFAULT);
 78:   EPSSetFromOptions(eps);
 79:   EPSGetST(eps,&st);
 80:   STGetKSP(st,&ksp);
 81:   KSPGetPC(ksp,&pc);
 82:   PCSetType(pc,PCSHELL);
 83:   PCShellSetApply(pc,PCApply_User);
 84:   STPrecondSetMatForPC(st,A);

 86:   /* set random initial vector */
 87:   MatCreateVecs(A,&v0,NULL);
 88:   PetscRandomCreate(PETSC_COMM_WORLD,&rand);
 89:   PetscRandomSetFromOptions(rand);
 90:   PetscOptionsGetInt(NULL,NULL,"-seed",&seed,NULL);
 91:   PetscRandomSetSeed(rand,seed);
 92:   PetscRandomSeed(rand);
 93:   VecSetRandom(v0,rand);
 94:   EPSSetInitialSpace(eps,1,&v0);
 95:   /* call the solver */
 96:   EPSSolve(eps);

 98:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 99:                     Display solution and clean up
100:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
101:   EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
102:   EPSDestroy(&eps);
103:   MatDestroy(&A);
104:   VecDestroy(&v0);
105:   PetscRandomDestroy(&rand);
106:   SlepcFinalize();
107:   return ierr;
108: }