Actual source code: test41.c

slepc-3.22.1 2024-10-28
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 interface to external package EVSL.\n\n"
 12:   "This is based on ex3.c. The command line options are:\n"
 13:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 14:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";

 16: #include <slepceps.h>

 18: int main(int argc,char **argv)
 19: {
 20:   Mat            A;           /* matrix */
 21:   EPS            eps;         /* eigenproblem solver context */
 22:   PetscInt       N,n=20,m,Istart,Iend,II,i,j,nslice;
 23:   PetscReal      a,b,ra,rb;
 24:   PetscBool      flag;
 25:   EPSEVSLDamping damping;

 27:   PetscFunctionBeginUser;
 28:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));

 30:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 31:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
 32:   if (!flag) m=n;
 33:   N = n*m;
 34:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nStandard eigenproblem with EVSL, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));

 36:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 37:      Compute the matrices that define the eigensystem, Ax=kBx
 38:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 40:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 41:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
 42:   PetscCall(MatSetFromOptions(A));

 44:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 45:   for (II=Istart;II<Iend;II++) {
 46:     i = II/n; j = II-i*n;
 47:     if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
 48:     if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
 49:     if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
 50:     if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
 51:     PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
 52:   }

 54:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 55:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));

 57:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 58:                 Create the eigensolver and set various options
 59:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 61:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
 62:   PetscCall(EPSSetOperators(eps,A,NULL));
 63:   PetscCall(EPSSetProblemType(eps,EPS_HEP));
 64:   PetscCall(EPSSetType(eps,EPSEVSL));

 66:   /*
 67:      Set several options
 68:   */
 69:   PetscCall(EPSSetInterval(eps,1.3,1.44));
 70:   PetscCall(EPSEVSLSetRange(eps,0,8));
 71:   PetscCall(EPSEVSLSetSlices(eps,3));
 72:   PetscCall(EPSEVSLSetDamping(eps,EPS_EVSL_DAMPING_SIGMA));
 73:   PetscCall(EPSEVSLSetDOSParameters(eps,EPS_EVSL_DOS_KPM,50,450,PETSC_CURRENT,PETSC_CURRENT));
 74:   PetscCall(EPSEVSLSetPolParameters(eps,4000,0.85));
 75:   PetscCall(EPSSetFromOptions(eps));

 77:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 78:            Compute all eigenvalues in interval and display info
 79:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 81:   PetscCall(EPSSolve(eps));
 82:   PetscCall(EPSEVSLGetSlices(eps,&nslice));
 83:   PetscCall(EPSGetInterval(eps,&a,&b));
 84:   PetscCall(EPSEVSLGetRange(eps,&ra,&rb));
 85:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," EVSL: solving interval [%g,%g] with %" PetscInt_FMT " slices (spectral range [%g,%g])\n",a,b,nslice,ra,rb));
 86:   PetscCall(EPSEVSLGetDamping(eps,&damping));
 87:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," EVSL: damping type is %s\n",EPSEVSLDampings[damping]));

 89:   PetscCall(EPSView(eps,NULL));
 90:   PetscCall(EPSErrorView(eps,EPS_ERROR_ABSOLUTE,NULL));

 92:   PetscCall(EPSDestroy(&eps));
 93:   PetscCall(MatDestroy(&A));
 94:   PetscCall(SlepcFinalize());
 95:   return 0;
 96: }

 98: /*TEST

100:    build:
101:       requires: evsl

103:    test:
104:       requires: evsl
105:       args: -eps_evsl_damping none
106:       filter: grep -v ncv

108: TEST*/