Actual source code: test26.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[] = "Illustrates the PGNHEP problem type. "
 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:   Mat               A,B;
 22:   PetscBool         flg;
 23:   PetscReal         tol=1000*PETSC_MACHINE_EPSILON;
 24:   char              filename[PETSC_MAX_PATH_LEN];
 25:   PetscViewer       viewer;

 27:   PetscFunctionBeginUser;
 28:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
 29:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nPGNHEP problem loaded from file\n\n"));

 31:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 32:         Load the matrices that define the eigensystem, Ax=kBx
 33:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

 61:   /* This example is intended for a matrix pair (A,B) where B is symmetric positive definite;
 62:      If we load matrices bfw62a/bfw62b, scale both of them because bfw62b is negative definite */
 63:   PetscCall(PetscStrendswith(filename,"bfw62b.petsc",&flg));
 64:   if (flg) {
 65:     PetscCall(MatScale(A,-1.0));
 66:     PetscCall(MatScale(B,-1.0));
 67:   }

 69:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 70:                 Create the eigensolver and solve the problem
 71:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 73:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
 74:   PetscCall(EPSSetOperators(eps,A,B));
 75:   PetscCall(EPSSetProblemType(eps,EPS_PGNHEP));
 76:   PetscCall(EPSSetTolerances(eps,tol,PETSC_CURRENT));
 77:   PetscCall(EPSSetFromOptions(eps));
 78:   PetscCall(EPSSolve(eps));

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

 84:   PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
 85:   PetscCall(EPSDestroy(&eps));
 86:   PetscCall(MatDestroy(&A));
 87:   PetscCall(MatDestroy(&B));
 88:   PetscCall(SlepcFinalize());
 89:   return 0;
 90: }

 92: /*TEST

 94:    testset:
 95:       args: -f1 ${SLEPC_DIR}/share/slepc/datafiles/matrices/bfw62a.petsc -f2 ${SLEPC_DIR}/share/slepc/datafiles/matrices/bfw62b.petsc -eps_largest_real -eps_nev 4
 96:       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
 97:       output_file: output/test26_1.out
 98:       test:
 99:          args: -eps_true_residual {{0 1}}
100:          suffix: 1
101:       test:
102:          args: -eps_type arpack
103:          suffix: 1_arpack
104:          requires: arpack

106:    testset:
107:       args: -f1 ${DATAFILESPATH}/matrices/complex/mhd1280a.petsc -f2 ${DATAFILESPATH}/matrices/complex/mhd1280b.petsc -eps_smallest_real -eps_nev 4
108:       requires: double complex datafilespath !defined(PETSC_USE_64BIT_INDICES)
109:       output_file: output/test26_2.out
110:       test:
111:          suffix: 2
112:       test:
113:          args: -eps_type arpack
114:          suffix: 2_arpack
115:          requires: arpack

117: TEST*/