Actual source code: test5.c
slepc-3.20.2 2024-03-15
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[] = "Test PEP view and monitor functionality.\n\n";
13: #include <slepcpep.h>
15: int main(int argc,char **argv)
16: {
17: Mat A[3];
18: PEP pep;
19: Vec xr,xi;
20: PetscScalar kr,ki;
21: PetscInt n=6,Istart,Iend,i,nconv,its;
22: PetscReal errest;
23: PetscBool checkfile;
24: char filename[PETSC_MAX_PATH_LEN];
25: PetscViewer viewer;
27: PetscFunctionBeginUser;
28: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
29: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
30: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nPEP of diagonal problem, n=%" PetscInt_FMT "\n\n",n));
32: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
33: Generate the matrices
34: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
35: PetscCall(MatCreate(PETSC_COMM_WORLD,&A[0]));
36: PetscCall(MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,n,n));
37: PetscCall(MatSetFromOptions(A[0]));
38: PetscCall(MatSetUp(A[0]));
39: PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend));
40: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(A[0],i,i,i+1,INSERT_VALUES));
41: PetscCall(MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY));
42: PetscCall(MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY));
44: PetscCall(MatCreate(PETSC_COMM_WORLD,&A[1]));
45: PetscCall(MatSetSizes(A[1],PETSC_DECIDE,PETSC_DECIDE,n,n));
46: PetscCall(MatSetFromOptions(A[1]));
47: PetscCall(MatSetUp(A[1]));
48: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(A[1],i,i,-1.5,INSERT_VALUES));
49: PetscCall(MatAssemblyBegin(A[1],MAT_FINAL_ASSEMBLY));
50: PetscCall(MatAssemblyEnd(A[1],MAT_FINAL_ASSEMBLY));
52: PetscCall(MatCreate(PETSC_COMM_WORLD,&A[2]));
53: PetscCall(MatSetSizes(A[2],PETSC_DECIDE,PETSC_DECIDE,n,n));
54: PetscCall(MatSetFromOptions(A[2]));
55: PetscCall(MatSetUp(A[2]));
56: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(A[2],i,i,-1.0/(i+1),INSERT_VALUES));
57: PetscCall(MatAssemblyBegin(A[2],MAT_FINAL_ASSEMBLY));
58: PetscCall(MatAssemblyEnd(A[2],MAT_FINAL_ASSEMBLY));
60: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61: Create the PEP solver
62: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63: PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
64: PetscCall(PetscObjectSetName((PetscObject)pep,"pep"));
65: PetscCall(PEPSetOperators(pep,3,A));
66: PetscCall(PEPSetFromOptions(pep));
68: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69: Solve the eigensystem and display solution
70: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71: PetscCall(PEPSolve(pep));
72: PetscCall(PEPGetConverged(pep,&nconv));
73: PetscCall(PEPGetIterationNumber(pep,&its));
74: PetscCall(PetscPrintf(PETSC_COMM_WORLD," %" PetscInt_FMT " converged eigenpairs after %" PetscInt_FMT " iterations\n",nconv,its));
75: if (nconv>0) {
76: PetscCall(MatCreateVecs(A[0],&xr,&xi));
77: PetscCall(PEPGetEigenpair(pep,0,&kr,&ki,xr,xi));
78: PetscCall(VecDestroy(&xr));
79: PetscCall(VecDestroy(&xi));
80: PetscCall(PEPGetErrorEstimate(pep,0,&errest));
81: }
82: PetscCall(PEPErrorView(pep,PEP_ERROR_RELATIVE,NULL));
84: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85: Check file containing the eigenvalues
86: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87: PetscCall(PetscOptionsGetString(NULL,NULL,"-checkfile",filename,sizeof(filename),&checkfile));
88: if (checkfile) {
89: #if defined(PETSC_HAVE_COMPLEX)
90: PetscComplex *eigs,eval;
91: PetscCall(PetscMalloc1(nconv,&eigs));
92: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer));
93: PetscCall(PetscViewerBinaryRead(viewer,eigs,nconv,NULL,PETSC_COMPLEX));
94: PetscCall(PetscViewerDestroy(&viewer));
95: for (i=0;i<nconv;i++) {
96: PetscCall(PEPGetEigenpair(pep,i,&kr,&ki,NULL,NULL));
97: #if defined(PETSC_USE_COMPLEX)
98: eval = kr;
99: #else
100: eval = PetscCMPLX(kr,ki);
101: #endif
102: PetscCheck(eval==eigs[i],PETSC_COMM_WORLD,PETSC_ERR_FILE_UNEXPECTED,"Eigenvalues in the file do not match");
103: }
104: PetscCall(PetscFree(eigs));
105: #else
106: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"The -checkfile option requires C99 complex numbers");
107: #endif
108: }
110: PetscCall(PEPDestroy(&pep));
111: PetscCall(MatDestroy(&A[0]));
112: PetscCall(MatDestroy(&A[1]));
113: PetscCall(MatDestroy(&A[2]));
114: PetscCall(SlepcFinalize());
115: return 0;
116: }
118: /*TEST
120: test:
121: suffix: 1
122: args: -pep_error_backward ::ascii_info_detail -pep_largest_real -pep_view_values -pep_monitor_conv -pep_error_absolute ::ascii_matlab -pep_monitor_all -pep_converged_reason -pep_view
123: requires: !single
124: filter: grep -v "tolerance" | grep -v "problem type" | sed -e "s/[+-]0\.0*i//g" -e "s/\([0-9]\.[5]*\)[+-][0-9]\.[0-9]*e-[0-9]*i/\\1/g" -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g"
126: test:
127: suffix: 2
128: args: -n 12 -pep_largest_real -pep_monitor -pep_view_values ::ascii_matlab
129: requires: double
130: filter: sed -e "s/[+-][0-9]\.[0-9]*e-[0-9]*i//" -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g" -e "s/5\.\([49]\)999999[0-9]*e+00/5.\\1999999999999999e+00/"
132: test:
133: suffix: 3
134: args: -pep_nev 4 -pep_view_values binary:myvalues.bin -checkfile myvalues.bin
135: requires: double c99_complex
137: test:
138: suffix: 4
139: args: -pep_nev 4 -pep_ncv 10 -pep_refine -pep_conv_norm -pep_extract none -pep_scale scalar -pep_view -pep_monitor -pep_error_relative ::ascii_info_detail
140: requires: double !complex
141: filter: grep -v "tolerance" | sed -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g"
143: test:
144: suffix: 5
145: args: -n 12 -pep_largest_real -pep_monitor draw::draw_lg -pep_monitor_all draw::draw_lg -pep_view_values draw -draw_save myeigen.ppm -draw_virtual
146: requires: x double
148: TEST*/