Actual source code: test41.c

slepc-3.15.1 2021-05-28
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 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;

 28:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

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

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

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

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

 55:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 56:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

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

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

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

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

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

 90:   EPSView(eps,NULL);
 91:   EPSErrorView(eps,EPS_ERROR_ABSOLUTE,NULL);

 93:   EPSDestroy(&eps);
 94:   MatDestroy(&A);
 95:   SlepcFinalize();
 96:   return ierr;
 97: }

 99: /*TEST

101:    build:
102:       requires: evsl

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

109: TEST*/