Actual source code: test26.c

slepc-3.17.0 2022-03-31
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:   SlepcInitialize(&argc,&argv,(char*)0,help);
 28:   PetscPrintf(PETSC_COMM_WORLD,"\nPGNHEP problem loaded from file\n\n");

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

 34:   PetscOptionsGetString(NULL,NULL,"-f1",filename,sizeof(filename),&flg);

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

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

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

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

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

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

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

 91: /*TEST

 93:    testset:
 94:       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
 95:       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
 96:       output_file: output/test26_1.out
 97:       test:
 98:          args: -eps_true_residual {{0 1}}
 99:          suffix: 1
100:       test:
101:          args: -eps_type arpack
102:          suffix: 1_arpack
103:          requires: arpack

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

116: TEST*/