Actual source code: test4.c

slepc-3.19.0 2023-03-31
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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 the solution of a HEP without calling EPSSetFromOptions (based on ex1.c).\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = number of grid subdivisions = matrix dimension.\n"
 14:   "  -type <eps_type> = eps type to test.\n\n";

 16: #include <slepceps.h>

 18: int main(int argc,char **argv)
 19: {
 20:   Mat            A;           /* problem matrix */
 21:   EPS            eps;         /* eigenproblem solver context */
 22:   PetscReal      tol=1000*PETSC_MACHINE_EPSILON;
 23:   PetscInt       n=30,i,Istart,Iend,nev;
 24:   PetscBool      flg,gd2;
 25:   char           epstype[30] = "krylovschur";

 27:   PetscFunctionBeginUser;
 28:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));

 30:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 31:   PetscCall(PetscOptionsGetString(NULL,NULL,"-type",epstype,sizeof(epstype),NULL));
 32:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian Eigenproblem, n=%" PetscInt_FMT "\n\n",n));

 34:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 35:      Compute the operator matrix that defines the eigensystem, Ax=kx
 36:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 38:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 39:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
 40:   PetscCall(MatSetFromOptions(A));
 41:   PetscCall(MatSetUp(A));

 43:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 44:   for (i=Istart;i<Iend;i++) {
 45:     if (i>0) PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
 46:     if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
 47:     PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES));
 48:   }
 49:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 50:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));

 52:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 53:                 Create the eigensolver and set various options
 54:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 55:   /*
 56:      Create eigensolver context
 57:   */
 58:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));

 60:   /*
 61:      Set operators. In this case, it is a standard eigenvalue problem
 62:   */
 63:   PetscCall(EPSSetOperators(eps,A,NULL));
 64:   PetscCall(EPSSetProblemType(eps,EPS_HEP));
 65:   PetscCall(EPSSetDimensions(eps,4,PETSC_DEFAULT,PETSC_DEFAULT));
 66:   PetscCall(EPSSetTolerances(eps,tol,PETSC_DEFAULT));

 68:   /*
 69:      Set solver parameters at runtime
 70:   */
 71:   PetscCall(PetscStrcmp(epstype,"gd2",&flg));
 72:   if (flg) {
 73:     PetscCall(EPSSetType(eps,EPSGD));
 74:     PetscCall(EPSGDSetDoubleExpansion(eps,PETSC_TRUE));
 75:     PetscCall(EPSGDGetDoubleExpansion(eps,&gd2));  /* not used */
 76:   } else PetscCall(EPSSetType(eps,epstype));
 77:   PetscCall(PetscStrcmp(epstype,"jd",&flg));
 78:   if (flg) {
 79:     PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
 80:     PetscCall(EPSSetTarget(eps,4.0));
 81:   }
 82:   PetscCall(PetscStrcmp(epstype,"lanczos",&flg));
 83:   if (flg) PetscCall(EPSLanczosSetReorthog(eps,EPS_LANCZOS_REORTHOG_LOCAL));
 84:   PetscCall(PetscObjectTypeCompareAny((PetscObject)eps,&flg,EPSRQCG,EPSLOBPCG,""));
 85:   if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));

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

 91:   PetscCall(EPSSolve(eps));
 92:   PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
 93:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));

 95:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 96:                     Display solution and clean up
 97:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 99:   PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
100:   PetscCall(EPSDestroy(&eps));
101:   PetscCall(MatDestroy(&A));
102:   PetscCall(SlepcFinalize());
103:   return 0;
104: }

106: /*TEST

108:    testset:
109:       output_file: output/test4_1.out
110:       test:
111:          suffix: 1
112:          args: -type {{krylovschur subspace arnoldi lanczos gd jd gd2 lapack}}
113:       test:
114:          suffix: 1_arpack
115:          args: -type arpack
116:          requires: arpack
117:       test:
118:          suffix: 1_primme
119:          args: -type primme
120:          requires: primme
121:       test:
122:          suffix: 1_trlan
123:          args: -type trlan
124:          requires: trlan

126:    test:
127:       suffix: 2
128:       args: -type {{rqcg lobpcg}}

130: TEST*/