Actual source code: test25.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[] = "Solves a GNHEP problem with contour integral. "
 12:   "Based on ex7.\n"
 13:   "The command line options are:\n"
 14:   "  -f1 <filename> -f2 <filename>, PETSc binary files containing A and B.\n\n";

 16: #include <slepceps.h>

 18: int main(int argc,char **argv)
 19: {
 20:   EPS               eps;
 21:   RG                rg;
 22:   Mat               A,B;
 23:   PetscBool         flg;
 24:   EPSCISSExtraction extr;
 25:   EPSCISSQuadRule   quad;
 26:   char              filename[PETSC_MAX_PATH_LEN];
 27:   PetscViewer       viewer;

 29:   PetscFunctionBeginUser;
 30:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
 31:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nGNHEP problem with contour integral\n\n"));

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

 37:   PetscCall(PetscOptionsGetString(NULL,NULL,"-f1",filename,sizeof(filename),&flg));
 38:   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Must indicate a file name for matrix A with the -f1 option");

 40: #if defined(PETSC_USE_COMPLEX)
 41:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Reading COMPLEX matrices from binary files...\n"));
 42: #else
 43:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Reading REAL matrices from binary files...\n"));
 44: #endif
 45:   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer));
 46:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 47:   PetscCall(MatSetFromOptions(A));
 48:   PetscCall(MatLoad(A,viewer));
 49:   PetscCall(PetscViewerDestroy(&viewer));

 51:   PetscCall(PetscOptionsGetString(NULL,NULL,"-f2",filename,sizeof(filename),&flg));
 52:   if (flg) {
 53:     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer));
 54:     PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
 55:     PetscCall(MatSetFromOptions(B));
 56:     PetscCall(MatLoad(B,viewer));
 57:     PetscCall(PetscViewerDestroy(&viewer));
 58:   } else {
 59:     PetscCall(PetscPrintf(PETSC_COMM_WORLD," Matrix B was not provided, setting B=I\n\n"));
 60:     B = NULL;
 61:   }

 63:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 64:                 Create the eigensolver and solve the problem
 65:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 67:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
 68:   PetscCall(EPSSetOperators(eps,A,B));
 69:   PetscCall(EPSSetProblemType(eps,EPS_GNHEP));
 70:   PetscCall(EPSSetTolerances(eps,1e-9,PETSC_CURRENT));

 72:   /* set CISS solver with various options */
 73:   PetscCall(EPSSetType(eps,EPSCISS));
 74:   PetscCall(EPSCISSSetExtraction(eps,EPS_CISS_EXTRACTION_HANKEL));
 75:   PetscCall(EPSCISSSetQuadRule(eps,EPS_CISS_QUADRULE_CHEBYSHEV));
 76:   PetscCall(EPSCISSSetUseST(eps,PETSC_TRUE));
 77:   PetscCall(EPSGetRG(eps,&rg));
 78:   PetscCall(RGSetType(rg,RGINTERVAL));
 79:   PetscCall(RGIntervalSetEndpoints(rg,-3000.0,0.0,0.0,0.0));

 81:   PetscCall(EPSSetFromOptions(eps));

 83:   PetscCall(EPSSolve(eps));
 84:   PetscCall(PetscObjectTypeCompare((PetscObject)eps,EPSCISS,&flg));
 85:   if (flg) {
 86:     PetscCall(EPSCISSGetExtraction(eps,&extr));
 87:     PetscCall(EPSCISSGetQuadRule(eps,&quad));
 88:     PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solved with CISS using %s extraction with %s quadrature rule\n\n",EPSCISSExtractions[extr],EPSCISSQuadRules[quad]));
 89:   }

 91:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 92:                     Display solution and clean up
 93:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 95:   PetscCall(EPSErrorView(eps,EPS_ERROR_BACKWARD,NULL));
 96:   PetscCall(EPSDestroy(&eps));
 97:   PetscCall(MatDestroy(&A));
 98:   PetscCall(MatDestroy(&B));
 99:   PetscCall(SlepcFinalize());
100:   return 0;
101: }

103: /*TEST

105:    testset:
106:       args: -f1 ${SLEPC_DIR}/share/slepc/datafiles/matrices/bfw62a.petsc -f2 ${SLEPC_DIR}/share/slepc/datafiles/matrices/bfw62b.petsc
107:       output_file: output/test25_1.out
108:       test:
109:          suffix: 1
110:          requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
111:       test:
112:          suffix: 1_cuda
113:          args: -mat_type aijcusparse -st_pc_factor_mat_solver_type cusparse
114:          requires: cuda double !complex !defined(PETSC_USE_64BIT_INDICES)
115:       test:
116:          suffix: 1_hip
117:          args: -mat_type aijhipsparse -st_pc_factor_mat_solver_type hipsparse
118:          requires: hip double !complex !defined(PETSC_USE_64BIT_INDICES)

120:    testset:
121:       args: -f1 ${DATAFILESPATH}/matrices/complex/qc324.petsc -rg_type ellipse -rg_ellipse_center 1-0.09i -rg_ellipse_radius 0.15 -rg_ellipse_vscale 0.1
122:       output_file: output/test25_2.out
123:       test:
124:          suffix: 2
125:          requires: double complex datafilespath !defined(PETSC_USE_64BIT_INDICES)
126:       test:
127:          suffix: 2_cuda
128:          args: -mat_type aijcusparse -st_pc_factor_mat_solver_type cusparse
129:          requires: cuda double complex datafilespath !defined(PETSC_USE_64BIT_INDICES)
130:       test:
131:          suffix: 2_hip
132:          args: -mat_type aijhipsparse -st_pc_factor_mat_solver_type hipsparse
133:          requires: hip double complex datafilespath !defined(PETSC_USE_64BIT_INDICES)

135: TEST*/