Actual source code: test26.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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;
 26:   PetscErrorCode    ierr;

 28:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 29:   PetscPrintf(PETSC_COMM_WORLD,"\nPGNHEP problem loaded from file\n\n");

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

 35:   PetscOptionsGetString(NULL,NULL,"-f1",filename,PETSC_MAX_PATH_LEN,&flg);
 36:   if (!flg) SETERRQ(PETSC_COMM_WORLD,1,"Must indicate a file name for matrix A with the -f1 option");

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

 49:   PetscOptionsGetString(NULL,NULL,"-f2",filename,PETSC_MAX_PATH_LEN,&flg);
 50:   if (flg) {
 51:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
 52:     MatCreate(PETSC_COMM_WORLD,&B);
 53:     MatSetFromOptions(B);
 54:     MatLoad(B,viewer);
 55:     PetscViewerDestroy(&viewer);
 56:   } else {
 57:     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:      We will load matrices bfw62a/bfw62b, and scale both of them because bfw62b is negative definite */
 63:   MatScale(A,-1.0);
 64:   MatScale(B,-1.0);

 66:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 67:                 Create the eigensolver and solve the problem
 68:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

 81:   EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
 82:   EPSDestroy(&eps);
 83:   MatDestroy(&A);
 84:   MatDestroy(&B);
 85:   SlepcFinalize();
 86:   return ierr;
 87: }

 89: /*TEST

 91:    test:
 92:       suffix: 1
 93:       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 -eps_true_residual {{0 1}}
 94:       requires: double !complex !define(PETSC_USE_64BIT_INDICES)
 95:       output_file: output/test26_1.out

 97: TEST*/