Actual source code: ex17.c
slepc-3.22.1 2024-10-28
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: PetscFunctionBeginUser;
30: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
32: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
33: Load the matrices that define the polynomial eigenproblem
34: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
36: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nPolynomial eigenproblem stored in file.\n\n"));
37: #if defined(PETSC_USE_COMPLEX)
38: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Reading COMPLEX matrices from binary files...\n"));
39: #else
40: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Reading REAL matrices from binary files...\n"));
41: #endif
42: PetscCall(PetscOptionsGetStringArray(NULL,NULL,"-A",filenames,&nmat,&flg));
43: PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Must indicate a comma-separated list of file names with the -A option");
44: for (i=0;i<nmat;i++) {
45: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filenames[i],FILE_MODE_READ,&viewer));
46: PetscCall(MatCreate(PETSC_COMM_WORLD,&A[i]));
47: PetscCall(MatSetFromOptions(A[i]));
48: PetscCall(MatLoad(A[i],viewer));
49: PetscCall(PetscViewerDestroy(&viewer));
50: }
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Create the eigensolver and set various options
53: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: /*
56: Create eigensolver context
57: */
58: PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
60: /*
61: Set matrices
62: */
63: PetscCall(PEPSetOperators(pep,nmat,A));
64: /*
65: Set solver parameters at runtime
66: */
67: PetscCall(PEPSetFromOptions(pep));
69: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
70: Solve the eigensystem
71: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73: PetscCall(PEPSolve(pep));
74: PetscCall(PEPGetIterationNumber(pep,&its));
75: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %" PetscInt_FMT "\n",its));
77: /*
78: Optional: Get some information from the solver and display it
79: */
80: PetscCall(PEPGetDimensions(pep,&nev,NULL,NULL));
81: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
82: PetscCall(PEPGetTolerances(pep,&tol,&maxit));
83: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%" PetscInt_FMT "\n",(double)tol,maxit));
85: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86: Display solution and clean up
87: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
89: /* show detailed info unless -terse option is given by user */
90: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
91: if (terse) PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
92: else {
93: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
94: PetscCall(PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD));
95: PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD));
96: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
97: }
98: PetscCall(PEPDestroy(&pep));
99: for (i=0;i<nmat;i++) {
100: PetscCall(MatDestroy(&A[i]));
101: PetscCall(PetscFree(filenames[i]));
102: }
103: PetscCall(SlepcFinalize());
104: return 0;
105: }
107: /*TEST
109: test:
110: suffix: 1
111: 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
112: requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
114: TEST*/