Actual source code: ex17.c

slepc-3.17.1 2022-04-11
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 polynomial eigenproblem P(l)x = 0 with matrices loaded from a file.\n\n"
 12:   "The command line options are:\n"
 13:   "-A <filename1,filename2, ...> , where <filename1,.. > = matrices A0 ... files in PETSc binary form.\n\n";

 15: #include <slepcpep.h>

 17: #define MAX_MATRICES 40

 19: int main(int argc,char **argv)
 20: {
 21:   Mat            A[MAX_MATRICES]; /* problem matrices */
 22:   PEP            pep;             /* polynomial eigenproblem solver context */
 23:   PetscReal      tol;
 24:   PetscInt       nev,maxit,its,nmat=MAX_MATRICES,i;
 25:   char*          filenames[MAX_MATRICES];
 26:   PetscViewer    viewer;
 27:   PetscBool      flg,terse;

 29:   SlepcInitialize(&argc,&argv,(char*)0,help);

 31:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 32:         Load the matrices that define the polynomial eigenproblem
 33:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 35:   PetscPrintf(PETSC_COMM_WORLD,"\nPolynomial eigenproblem stored in file.\n\n");
 36: #if defined(PETSC_USE_COMPLEX)
 37:   PetscPrintf(PETSC_COMM_WORLD," Reading COMPLEX matrices from binary files...\n");
 38: #else
 39:   PetscPrintf(PETSC_COMM_WORLD," Reading REAL matrices from binary files...\n");
 40: #endif
 41:   PetscOptionsGetStringArray(NULL,NULL,"-A",filenames,&nmat,&flg);
 43:   for (i=0;i<nmat;i++) {
 44:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,filenames[i],FILE_MODE_READ,&viewer);
 45:     MatCreate(PETSC_COMM_WORLD,&A[i]);
 46:     MatSetFromOptions(A[i]);
 47:     MatLoad(A[i],viewer);
 48:     PetscViewerDestroy(&viewer);
 49:   }
 50:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 51:                 Create the eigensolver and set various options
 52:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 54:   /*
 55:      Create eigensolver context
 56:   */
 57:   PEPCreate(PETSC_COMM_WORLD,&pep);

 59:   /*
 60:      Set matrices
 61:   */
 62:   PEPSetOperators(pep,nmat,A);
 63:   /*
 64:      Set solver parameters at runtime
 65:   */
 66:   PEPSetFromOptions(pep);

 68:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 69:                       Solve the eigensystem
 70:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 72:   PEPSolve(pep);
 73:   PEPGetIterationNumber(pep,&its);
 74:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %" PetscInt_FMT "\n",its);

 76:   /*
 77:      Optional: Get some information from the solver and display it
 78:   */
 79:   PEPGetDimensions(pep,&nev,NULL,NULL);
 80:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
 81:   PEPGetTolerances(pep,&tol,&maxit);
 82:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%" PetscInt_FMT "\n",(double)tol,maxit);

 84:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 85:                     Display solution and clean up
 86:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 88:   /* show detailed info unless -terse option is given by user */
 89:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
 90:   if (terse) PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
 91:   else {
 92:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
 93:     PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
 94:     PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
 95:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
 96:   }
 97:   PEPDestroy(&pep);
 98:   for (i=0;i<nmat;i++) {
 99:     MatDestroy(&A[i]);
100:     PetscFree(filenames[i]);
101:   }
102:   SlepcFinalize();
103:   return 0;
104: }

106: /*TEST

108:    test:
109:       suffix: 1
110:       args: -A ${SLEPC_DIR}/share/slepc/datafiles/matrices/speaker107k.petsc,${SLEPC_DIR}/share/slepc/datafiles/matrices/speaker107c.petsc,${SLEPC_DIR}/share/slepc/datafiles/matrices/speaker107m.petsc -pep_type {{toar qarnoldi linear}} -pep_nev 4 -pep_ncv 20 -pep_scale scalar -terse
111:       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)

113: TEST*/