Actual source code: test8.c
slepc-3.17.2 2022-08-09
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 NEP view and monitor functionality.\n\n";
13: #include <slepcnep.h>
15: int main(int argc,char **argv)
16: {
17: Mat A[3];
18: FN f[3];
19: NEP nep;
20: Vec xr,xi;
21: PetscScalar kr,ki,coeffs[3];
22: PetscInt n=6,i,Istart,Iend,nconv,its;
23: PetscReal errest;
24: PetscBool checkfile;
25: char filename[PETSC_MAX_PATH_LEN];
26: PetscViewer viewer;
28: SlepcInitialize(&argc,&argv,(char*)0,help);
29: PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Nonlinear Eigenproblem, n=%" PetscInt_FMT "\n\n",n);
31: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
32: Generate the matrices
33: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
35: MatCreate(PETSC_COMM_WORLD,&A[0]);
36: MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,n,n);
37: MatSetFromOptions(A[0]);
38: MatSetUp(A[0]);
39: MatGetOwnershipRange(A[0],&Istart,&Iend);
40: for (i=Istart;i<Iend;i++) MatSetValue(A[0],i,i,i+1,INSERT_VALUES);
41: MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY);
42: MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY);
44: MatCreate(PETSC_COMM_WORLD,&A[1]);
45: MatSetSizes(A[1],PETSC_DECIDE,PETSC_DECIDE,n,n);
46: MatSetFromOptions(A[1]);
47: MatSetUp(A[1]);
48: for (i=Istart;i<Iend;i++) MatSetValue(A[1],i,i,-1.5,INSERT_VALUES);
49: MatAssemblyBegin(A[1],MAT_FINAL_ASSEMBLY);
50: MatAssemblyEnd(A[1],MAT_FINAL_ASSEMBLY);
52: MatCreate(PETSC_COMM_WORLD,&A[2]);
53: MatSetSizes(A[2],PETSC_DECIDE,PETSC_DECIDE,n,n);
54: MatSetFromOptions(A[2]);
55: MatSetUp(A[2]);
56: for (i=Istart;i<Iend;i++) MatSetValue(A[2],i,i,-1.0/(i+1),INSERT_VALUES);
57: MatAssemblyBegin(A[2],MAT_FINAL_ASSEMBLY);
58: MatAssemblyEnd(A[2],MAT_FINAL_ASSEMBLY);
60: /*
61: Functions: f0=1.0, f1=lambda, f2=lambda^2
62: */
63: FNCreate(PETSC_COMM_WORLD,&f[0]);
64: FNSetType(f[0],FNRATIONAL);
65: coeffs[0] = 1.0;
66: FNRationalSetNumerator(f[0],1,coeffs);
68: FNCreate(PETSC_COMM_WORLD,&f[1]);
69: FNSetType(f[1],FNRATIONAL);
70: coeffs[0] = 1.0; coeffs[1] = 0.0;
71: FNRationalSetNumerator(f[1],2,coeffs);
73: FNCreate(PETSC_COMM_WORLD,&f[2]);
74: FNSetType(f[2],FNRATIONAL);
75: coeffs[0] = 1.0; coeffs[1] = 0.0; coeffs[2] = 0.0;
76: FNRationalSetNumerator(f[2],3,coeffs);
78: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79: Create the NEP solver
80: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
81: NEPCreate(PETSC_COMM_WORLD,&nep);
82: PetscObjectSetName((PetscObject)nep,"nep");
83: NEPSetSplitOperator(nep,3,A,f,SAME_NONZERO_PATTERN);
84: NEPSetTarget(nep,1.1);
85: NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE);
86: NEPSetFromOptions(nep);
88: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89: Solve the eigensystem and display solution
90: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91: NEPSolve(nep);
92: NEPGetConverged(nep,&nconv);
93: NEPGetIterationNumber(nep,&its);
94: PetscPrintf(PETSC_COMM_WORLD," %" PetscInt_FMT " converged eigenpairs after %" PetscInt_FMT " iterations\n",nconv,its);
95: if (nconv>0) {
96: MatCreateVecs(A[0],&xr,&xi);
97: NEPGetEigenpair(nep,0,&kr,&ki,xr,xi);
98: VecDestroy(&xr);
99: VecDestroy(&xi);
100: NEPGetErrorEstimate(nep,0,&errest);
101: }
102: NEPErrorView(nep,NEP_ERROR_BACKWARD,NULL);
104: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105: Check file containing the eigenvalues
106: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107: PetscOptionsGetString(NULL,NULL,"-checkfile",filename,sizeof(filename),&checkfile);
108: if (checkfile) {
109: #if defined(PETSC_HAVE_COMPLEX)
110: PetscComplex *eigs,eval;
111: PetscMalloc1(nconv,&eigs);
112: PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
113: PetscViewerBinaryRead(viewer,eigs,nconv,NULL,PETSC_COMPLEX);
114: PetscViewerDestroy(&viewer);
115: for (i=0;i<nconv;i++) {
116: NEPGetEigenpair(nep,i,&kr,&ki,NULL,NULL);
117: #if defined(PETSC_USE_COMPLEX)
118: eval = kr;
119: #else
120: eval = PetscCMPLX(kr,ki);
121: #endif
123: }
124: PetscFree(eigs);
125: #else
126: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"The -checkfile option requires C99 complex numbers");
127: #endif
128: }
130: NEPDestroy(&nep);
131: MatDestroy(&A[0]);
132: MatDestroy(&A[1]);
133: MatDestroy(&A[2]);
134: FNDestroy(&f[0]);
135: FNDestroy(&f[1]);
136: FNDestroy(&f[2]);
137: SlepcFinalize();
138: return 0;
139: }
141: /*TEST
143: test:
144: suffix: 1
145: args: -nep_type slp -nep_target -.5 -nep_error_backward ::ascii_info_detail -nep_view_values -nep_error_absolute ::ascii_matlab -nep_monitor_all -nep_converged_reason -nep_view
146: filter: grep -v "tolerance" | grep -v "problem type" | sed -e "s/[+-]0\.0*i//g" -e "s/+0i//" -e "s/[+-][0-9]\.[0-9]*e-[0-9]*i//g" -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
147: requires: double
149: test:
150: suffix: 2
151: args: -nep_type rii -nep_target -.5 -nep_rii_hermitian -nep_monitor -nep_view_values ::ascii_matlab
152: filter: sed -e "s/[+-][0-9]\.[0-9]*e-[0-9]*i//" -e "s/([0-9]\.[0-9]*e[+-]\([0-9]*\))/(removed)/g"
153: requires: double
155: test:
156: suffix: 3
157: args: -nep_type slp -nep_nev 4 -nep_view_values binary:myvalues.bin -checkfile myvalues.bin -nep_error_relative ::ascii_matlab
158: filter: sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
159: requires: double c99_complex
161: test:
162: suffix: 4
163: args: -nep_type slp -nep_nev 4 -nep_monitor draw::draw_lg -nep_monitor_all draw::draw_lg -nep_view_values draw -draw_save myeigen.ppm -draw_virtual
164: requires: x double
166: TEST*/