Actual source code: test18.c

slepc-3.16.2 2022-02-01
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[] = "Symmetric-indefinite eigenproblem.\n\n"
 12:   "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,B;             /* problem matrices */
 21:   EPS            eps;             /* eigenproblem solver context */
 22:   PetscInt       N,n=10,m,Istart,Iend,II,nev,i,j;
 23:   PetscBool      flag,terse;

 26:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 27:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 28:   PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
 29:   if (!flag) m=n;
 30:   N = n*m;
 31:   PetscPrintf(PETSC_COMM_WORLD,"\nSymmetric-indefinite eigenproblem, N=%D\n\n",N);

 33:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 34:           Compute the matrices that define the eigensystem, Ax=kBx
 35:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 37:   MatCreate(PETSC_COMM_WORLD,&A);
 38:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 39:   MatSetFromOptions(A);
 40:   MatSetUp(A);

 42:   MatCreate(PETSC_COMM_WORLD,&B);
 43:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);
 44:   MatSetFromOptions(B);
 45:   MatSetUp(B);

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

 58:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 59:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 60:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 61:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 63:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 64:                 Create the eigensolver and set various options
 65:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 67:   EPSCreate(PETSC_COMM_WORLD,&eps);
 68:   EPSSetOperators(eps,A,B);
 69:   EPSSetProblemType(eps,EPS_GHIEP);
 70:   EPSSetFromOptions(eps);

 72:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 73:                       Solve the eigensystem
 74:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 76:   EPSSolve(eps);
 77:   EPSGetDimensions(eps,&nev,NULL,NULL);
 78:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

 80:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 81:                     Display solution and clean up
 82:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 84:   /* show detailed info unless -terse option is given by user */
 85:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
 86:   if (terse) {
 87:     EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
 88:   } else {
 89:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
 90:     EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
 91:     EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
 92:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
 93:   }

 95:   EPSDestroy(&eps);
 96:   MatDestroy(&A);
 97:   MatDestroy(&B);
 98:   SlepcFinalize();
 99:   return ierr;
100: }

102: /*TEST

104:    testset:
105:       args: -eps_nev 4 -eps_ncv 12 -terse -st_type sinvert -eps_krylovschur_restart .3
106:       requires: !single
107:       output_file: output/test18_1.out
108:       test:
109:          suffix: 1_ks
110:       test:
111:          suffix: 1_ks_gnhep
112:          args: -eps_gen_non_hermitian
113:       test:
114:          suffix: 2_cuda_ks
115:          args: -mat_type aijcusparse
116:          requires: cuda
117:       test:
118:          suffix: 2_cuda_ks_gnhep
119:          args: -eps_gen_non_hermitian -mat_type aijcusparse
120:          requires: cuda

122:    test:
123:       suffix: 2
124:       args: -n 10 -m 11 -eps_type {{gd jd}} -eps_target 0.2 -eps_harmonic -eps_nev 2 -eps_ncv 10 -terse
125:       requires: !single
126:       filter: sed -e "s/[+-]0\.0*i//g"

128: TEST*/