LCOV - code coverage report
Current view: top level - eps/tests - test41.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 50 50 100.0 %
Date: 2024-05-08 00:46:20 Functions: 1 1 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /*
       2             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       3             :    SLEPc - Scalable Library for Eigenvalue Problem Computations
       4             :    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
       5             : 
       6             :    This file is part of SLEPc.
       7             :    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
       8             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       9             : */
      10             : 
      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";
      15             : 
      16             : #include <slepceps.h>
      17             : 
      18           1 : int main(int argc,char **argv)
      19             : {
      20           1 :   Mat            A;           /* matrix */
      21           1 :   EPS            eps;         /* eigenproblem solver context */
      22           1 :   PetscInt       N,n=20,m,Istart,Iend,II,i,j,nslice;
      23           1 :   PetscReal      a,b,ra,rb;
      24           1 :   PetscBool      flag;
      25           1 :   EPSEVSLDamping damping;
      26             : 
      27           1 :   PetscFunctionBeginUser;
      28           1 :   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
      29             : 
      30           1 :   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
      31           1 :   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
      32           1 :   if (!flag) m=n;
      33           1 :   N = n*m;
      34           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nStandard eigenproblem with EVSL, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
      35             : 
      36             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
      37             :      Compute the matrices that define the eigensystem, Ax=kBx
      38             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
      39             : 
      40           1 :   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
      41           1 :   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
      42           1 :   PetscCall(MatSetFromOptions(A));
      43             : 
      44           1 :   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
      45         401 :   for (II=Istart;II<Iend;II++) {
      46         400 :     i = II/n; j = II-i*n;
      47         400 :     if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
      48         400 :     if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
      49         400 :     if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
      50         400 :     if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
      51         400 :     PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
      52             :   }
      53             : 
      54           1 :   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
      55           1 :   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
      56             : 
      57             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
      58             :                 Create the eigensolver and set various options
      59             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
      60             : 
      61           1 :   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
      62           1 :   PetscCall(EPSSetOperators(eps,A,NULL));
      63           1 :   PetscCall(EPSSetProblemType(eps,EPS_HEP));
      64           1 :   PetscCall(EPSSetType(eps,EPSEVSL));
      65             : 
      66             :   /*
      67             :      Set several options
      68             :   */
      69           1 :   PetscCall(EPSSetInterval(eps,1.3,1.44));
      70           1 :   PetscCall(EPSEVSLSetRange(eps,0,8));
      71           1 :   PetscCall(EPSEVSLSetSlices(eps,3));
      72           1 :   PetscCall(EPSEVSLSetDamping(eps,EPS_EVSL_DAMPING_SIGMA));
      73           1 :   PetscCall(EPSEVSLSetDOSParameters(eps,EPS_EVSL_DOS_KPM,50,450,PETSC_DEFAULT,PETSC_DEFAULT));
      74           1 :   PetscCall(EPSEVSLSetPolParameters(eps,4000,0.85));
      75           1 :   PetscCall(EPSSetFromOptions(eps));
      76             : 
      77             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
      78             :            Compute all eigenvalues in interval and display info
      79             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
      80             : 
      81           1 :   PetscCall(EPSSolve(eps));
      82           1 :   PetscCall(EPSEVSLGetSlices(eps,&nslice));
      83           1 :   PetscCall(EPSGetInterval(eps,&a,&b));
      84           1 :   PetscCall(EPSEVSLGetRange(eps,&ra,&rb));
      85           1 :   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           1 :   PetscCall(EPSEVSLGetDamping(eps,&damping));
      87           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD," EVSL: damping type is %s\n",EPSEVSLDampings[damping]));
      88             : 
      89           1 :   PetscCall(EPSView(eps,NULL));
      90           1 :   PetscCall(EPSErrorView(eps,EPS_ERROR_ABSOLUTE,NULL));
      91             : 
      92           1 :   PetscCall(EPSDestroy(&eps));
      93           1 :   PetscCall(MatDestroy(&A));
      94           1 :   PetscCall(SlepcFinalize());
      95             :   return 0;
      96             : }
      97             : 
      98             : /*TEST
      99             : 
     100             :    build:
     101             :       requires: evsl
     102             : 
     103             :    test:
     104             :       requires: evsl
     105             :       args: -eps_evsl_damping none
     106             :       filter: grep -v ncv
     107             : 
     108             : TEST*/

Generated by: LCOV version 1.14