Actual source code: test8.c
slepc-3.15.1 2021-05-28
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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: PetscComplex *eigs,eval;
23: PetscInt n=6,i,Istart,Iend,nconv,its;
24: PetscReal errest;
25: PetscBool checkfile;
26: char filename[PETSC_MAX_PATH_LEN];
27: PetscViewer viewer;
28: PetscErrorCode ierr;
30: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
31: PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Nonlinear Eigenproblem, n=%D\n\n",n);
33: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34: Generate the matrices
35: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
37: MatCreate(PETSC_COMM_WORLD,&A[0]);
38: MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,n,n);
39: MatSetFromOptions(A[0]);
40: MatSetUp(A[0]);
41: MatGetOwnershipRange(A[0],&Istart,&Iend);
42: for (i=Istart;i<Iend;i++) {
43: MatSetValue(A[0],i,i,i+1,INSERT_VALUES);
44: }
45: MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY);
46: MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY);
48: MatCreate(PETSC_COMM_WORLD,&A[1]);
49: MatSetSizes(A[1],PETSC_DECIDE,PETSC_DECIDE,n,n);
50: MatSetFromOptions(A[1]);
51: MatSetUp(A[1]);
52: for (i=Istart;i<Iend;i++) {
53: MatSetValue(A[1],i,i,-1.5,INSERT_VALUES);
54: }
55: MatAssemblyBegin(A[1],MAT_FINAL_ASSEMBLY);
56: MatAssemblyEnd(A[1],MAT_FINAL_ASSEMBLY);
58: MatCreate(PETSC_COMM_WORLD,&A[2]);
59: MatSetSizes(A[2],PETSC_DECIDE,PETSC_DECIDE,n,n);
60: MatSetFromOptions(A[2]);
61: MatSetUp(A[2]);
62: for (i=Istart;i<Iend;i++) {
63: MatSetValue(A[2],i,i,-1.0/(i+1),INSERT_VALUES);
64: }
65: MatAssemblyBegin(A[2],MAT_FINAL_ASSEMBLY);
66: MatAssemblyEnd(A[2],MAT_FINAL_ASSEMBLY);
68: /*
69: Functions: f0=1.0, f1=lambda, f2=lambda^2
70: */
71: FNCreate(PETSC_COMM_WORLD,&f[0]);
72: FNSetType(f[0],FNRATIONAL);
73: coeffs[0] = 1.0;
74: FNRationalSetNumerator(f[0],1,coeffs);
76: FNCreate(PETSC_COMM_WORLD,&f[1]);
77: FNSetType(f[1],FNRATIONAL);
78: coeffs[0] = 1.0; coeffs[1] = 0.0;
79: FNRationalSetNumerator(f[1],2,coeffs);
81: FNCreate(PETSC_COMM_WORLD,&f[2]);
82: FNSetType(f[2],FNRATIONAL);
83: coeffs[0] = 1.0; coeffs[1] = 0.0; coeffs[2] = 0.0;
84: FNRationalSetNumerator(f[2],3,coeffs);
86: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87: Create the NEP solver
88: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
89: NEPCreate(PETSC_COMM_WORLD,&nep);
90: PetscObjectSetName((PetscObject)nep,"nep");
91: NEPSetSplitOperator(nep,3,A,f,SAME_NONZERO_PATTERN);
92: NEPSetTarget(nep,1.1);
93: NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE);
94: NEPSetFromOptions(nep);
96: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: Solve the eigensystem and display solution
98: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99: NEPSolve(nep);
100: NEPGetConverged(nep,&nconv);
101: NEPGetIterationNumber(nep,&its);
102: PetscPrintf(PETSC_COMM_WORLD," %D converged eigenpairs after %D iterations\n",nconv,its);
103: if (nconv>0) {
104: MatCreateVecs(A[0],&xr,&xi);
105: NEPGetEigenpair(nep,0,&kr,&ki,xr,xi);
106: VecDestroy(&xr);
107: VecDestroy(&xi);
108: NEPGetErrorEstimate(nep,0,&errest);
109: }
110: NEPErrorView(nep,NEP_ERROR_BACKWARD,NULL);
112: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113: Check file containing the eigenvalues
114: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115: PetscOptionsGetString(NULL,NULL,"-checkfile",filename,sizeof(filename),&checkfile);
116: if (checkfile) {
117: PetscMalloc1(nconv,&eigs);
118: PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
119: PetscViewerBinaryRead(viewer,eigs,nconv,NULL,PETSC_COMPLEX);
120: PetscViewerDestroy(&viewer);
121: for (i=0;i<nconv;i++) {
122: NEPGetEigenpair(nep,i,&kr,&ki,NULL,NULL);
123: #if defined(PETSC_USE_COMPLEX)
124: eval = kr;
125: #else
126: eval = PetscCMPLX(kr,ki);
127: #endif
128: if (eval!=eigs[i]) SETERRQ(PETSC_COMM_WORLD,1,"Eigenvalues in the file do not match");
129: }
130: PetscFree(eigs);
131: }
133: NEPDestroy(&nep);
134: MatDestroy(&A[0]);
135: MatDestroy(&A[1]);
136: MatDestroy(&A[2]);
137: FNDestroy(&f[0]);
138: FNDestroy(&f[1]);
139: FNDestroy(&f[2]);
140: SlepcFinalize();
141: return ierr;
142: }
144: /*TEST
146: test:
147: suffix: 1
148: 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
149: 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"
150: requires: double
152: test:
153: suffix: 2
154: args: -nep_type rii -nep_target -.5 -nep_rii_hermitian -nep_monitor -nep_view_values ::ascii_matlab
155: filter: sed -e "s/[+-][0-9]\.[0-9]*e-[0-9]*i//" -e "s/([0-9]\.[0-9]*e[+-]\([0-9]*\))/(removed)/g"
156: requires: double
158: test:
159: suffix: 3
160: args: -nep_type slp -nep_nev 4 -nep_view_values binary:myvalues.bin -checkfile myvalues.bin -nep_error_relative ::ascii_matlab
161: filter: sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
162: requires: double
164: test:
165: suffix: 4
166: 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
167: requires: double
169: TEST*/