Actual source code: test18.c

slepc-3.22.2 2024-12-02
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[] = "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;

 25:   PetscFunctionBeginUser;
 26:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
 27:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 28:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
 29:   if (!flag) m=n;
 30:   N = n*m;
 31:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSymmetric-indefinite eigenproblem, N=%" PetscInt_FMT "\n\n",N));

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

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

 41:   PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
 42:   PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N));
 43:   PetscCall(MatSetFromOptions(B));

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

 56:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 57:   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
 58:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
 59:   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));

 61:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 62:                 Create the eigensolver and set various options
 63:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 65:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
 66:   PetscCall(EPSSetOperators(eps,A,B));
 67:   PetscCall(EPSSetProblemType(eps,EPS_GHIEP));
 68:   PetscCall(EPSSetFromOptions(eps));

 70:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 71:                       Solve the eigensystem
 72:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 74:   PetscCall(EPSSolve(eps));
 75:   PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
 76:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));

 78:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 79:                     Display solution and clean up
 80:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

 99: /*TEST

101:    testset:
102:       args: -eps_nev 4 -eps_ncv 12 -terse -st_type sinvert -eps_krylovschur_restart .3
103:       requires: !single
104:       output_file: output/test18_1.out
105:       test:
106:          suffix: 1_ks
107:       test:
108:          suffix: 1_ks_gnhep
109:          args: -eps_gen_non_hermitian
110:          requires: !__float128
111:       test:
112:          suffix: 2_cuda_ks
113:          args: -mat_type aijcusparse
114:          requires: cuda
115:       test:
116:          suffix: 2_cuda_ks_gnhep
117:          args: -eps_gen_non_hermitian -mat_type aijcusparse
118:          requires: cuda
119:       test:
120:          suffix: 2_hip_ks
121:          args: -mat_type aijhipsparse
122:          requires: hip
123:       test:
124:          suffix: 2_hip_ks_gnhep
125:          args: -eps_gen_non_hermitian -mat_type aijhipsparse
126:          requires: hip

128:    testset:
129:       args: -n 10 -m 11 -eps_target 0.2 -eps_harmonic -eps_nev 2 -terse
130:       filter: sed -e "s/[+-]0\.0*i//g"
131:       output_file: output/test18_2.out
132:       test:
133:          suffix: 2_gd
134:          args: -eps_type gd -eps_ncv 8
135:          requires: !single
136:       test:
137:          suffix: 2_jd
138:          args: -eps_type jd -st_ksp_type bcgs -eps_ncv 10
139:          requires: double

141: TEST*/