Actual source code: test30.c

slepc-3.16.3 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 11: static char help[] = "Test changing EPS type.\n\n";

 13: #include <slepceps.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Mat            A;           /* problem matrix */
 18:   EPS            eps;         /* eigenproblem solver context */
 19:   ST             st;
 20:   KSP            ksp;
 21:   PetscInt       n=20,i,Istart,Iend,nrest;

 24:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 25:   PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Eigenproblem, n=%D\n\n",n);

 27:   MatCreate(PETSC_COMM_WORLD,&A);
 28:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 29:   MatSetFromOptions(A);
 30:   MatSetUp(A);
 31:   MatGetOwnershipRange(A,&Istart,&Iend);
 32:   for (i=Istart;i<Iend;i++) {
 33:     MatSetValue(A,i,i,i+1,INSERT_VALUES);
 34:   }
 35:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 36:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 38:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 39:              Create eigensolver and view options
 40:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 41:   EPSCreate(PETSC_COMM_WORLD,&eps);
 42:   EPSSetOperators(eps,A,NULL);
 43:   EPSSetProblemType(eps,EPS_HEP);
 44:   EPSSetTarget(eps,4.8);
 45:   EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);
 46:   EPSSetDimensions(eps,4,PETSC_DEFAULT,PETSC_DEFAULT);

 48:   EPSSetType(eps,EPSRQCG);
 49:   EPSRQCGSetReset(eps,10);
 50:   EPSSetFromOptions(eps);
 51:   EPSRQCGGetReset(eps,&nrest);
 52:   PetscPrintf(PETSC_COMM_WORLD,"RQCG reset parameter = %D\n\n",nrest);
 53:   EPSView(eps,NULL);

 55:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 56:              Change eigensolver type and solve problem
 57:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 58:   EPSSetType(eps,EPSGD);
 59:   EPSGetST(eps,&st);
 60:   STGetKSP(st,&ksp);
 61:   KSPSetType(ksp,KSPPREONLY);
 62:   EPSSolve(eps);
 63:   EPSView(eps,NULL);

 65:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 66:                     Display solution and clean up
 67:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 68:   EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
 69:   EPSDestroy(&eps);
 70:   MatDestroy(&A);
 71:   SlepcFinalize();
 72:   return ierr;
 73: }

 75: /*TEST

 77:    test:
 78:       suffix: 1
 79:       args: -eps_harmonic -eps_conv_norm
 80:       filter: grep -v tolerance | sed -e "s/symmetric/hermitian/"

 82: TEST*/