Actual source code: test14.c
slepc-3.19.0 2023-03-31
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 EPS interface functions.\n\n";
13: #include <slepceps.h>
15: int main(int argc,char **argv)
16: {
17: Mat A,B; /* problem matrix */
18: EPS eps; /* eigenproblem solver context */
19: ST st;
20: KSP ksp;
21: DS ds;
22: PetscReal cut,tol;
23: PetscScalar target;
24: PetscInt n=20,i,its,nev,ncv,mpd,Istart,Iend;
25: PetscBool flg,pur,track;
26: EPSConvergedReason reason;
27: EPSType type;
28: EPSExtraction extr;
29: EPSBalance bal;
30: EPSWhich which;
31: EPSConv conv;
32: EPSStop stop;
33: EPSProblemType ptype;
34: PetscViewerAndFormat *vf;
36: PetscFunctionBeginUser;
37: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
38: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Eigenproblem, n=%" PetscInt_FMT "\n\n",n));
40: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
41: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
42: PetscCall(MatSetFromOptions(A));
43: PetscCall(MatSetUp(A));
44: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
45: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(A,i,i,i+1,INSERT_VALUES));
46: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
47: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
49: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50: Create eigensolver and test interface functions
51: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
52: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
53: PetscCall(EPSSetOperators(eps,A,NULL));
54: PetscCall(EPSGetOperators(eps,&B,NULL));
55: PetscCall(MatView(B,NULL));
57: PetscCall(EPSSetType(eps,EPSKRYLOVSCHUR));
58: PetscCall(EPSGetType(eps,&type));
59: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Type set to %s\n",type));
61: PetscCall(EPSGetProblemType(eps,&ptype));
62: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Problem type before changing = %d",(int)ptype));
63: PetscCall(EPSSetProblemType(eps,EPS_HEP));
64: PetscCall(EPSGetProblemType(eps,&ptype));
65: PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %d.",(int)ptype));
66: PetscCall(EPSIsGeneralized(eps,&flg));
67: if (flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD," generalized"));
68: PetscCall(EPSIsHermitian(eps,&flg));
69: if (flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD," hermitian"));
70: PetscCall(EPSIsPositive(eps,&flg));
71: if (flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD," positive"));
73: PetscCall(EPSGetExtraction(eps,&extr));
74: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n Extraction before changing = %d",(int)extr));
75: PetscCall(EPSSetExtraction(eps,EPS_HARMONIC));
76: PetscCall(EPSGetExtraction(eps,&extr));
77: PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)extr));
79: PetscCall(EPSSetBalance(eps,EPS_BALANCE_ONESIDE,8,1e-6));
80: PetscCall(EPSGetBalance(eps,&bal,&its,&cut));
81: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Balance: %s, its=%" PetscInt_FMT ", cutoff=%g\n",EPSBalanceTypes[bal],its,(double)cut));
83: PetscCall(EPSSetPurify(eps,PETSC_FALSE));
84: PetscCall(EPSGetPurify(eps,&pur));
85: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Eigenvector purification: %s\n",pur?"on":"off"));
86: PetscCall(EPSGetTrackAll(eps,&track));
88: PetscCall(EPSSetTarget(eps,4.8));
89: PetscCall(EPSGetTarget(eps,&target));
90: PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
91: PetscCall(EPSGetWhichEigenpairs(eps,&which));
92: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Which = %d, target = %g\n",(int)which,(double)PetscRealPart(target)));
94: PetscCall(EPSSetDimensions(eps,4,PETSC_DEFAULT,PETSC_DEFAULT));
95: PetscCall(EPSGetDimensions(eps,&nev,&ncv,&mpd));
96: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Dimensions: nev=%" PetscInt_FMT ", ncv=%" PetscInt_FMT ", mpd=%" PetscInt_FMT "\n",nev,ncv,mpd));
98: PetscCall(EPSSetTolerances(eps,2.2e-4,200));
99: PetscCall(EPSGetTolerances(eps,&tol,&its));
100: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Tolerance = %.5f, max_its = %" PetscInt_FMT "\n",(double)tol,its));
102: PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_ABS));
103: PetscCall(EPSGetConvergenceTest(eps,&conv));
104: PetscCall(EPSSetStoppingTest(eps,EPS_STOP_BASIC));
105: PetscCall(EPSGetStoppingTest(eps,&stop));
106: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Convergence test = %d, stopping test = %d\n",(int)conv,(int)stop));
108: PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf));
109: PetscCall(EPSMonitorSet(eps,(PetscErrorCode (*)(EPS,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))EPSMonitorFirst,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy));
110: PetscCall(EPSMonitorCancel(eps));
112: PetscCall(EPSGetST(eps,&st));
113: PetscCall(STGetKSP(st,&ksp));
114: PetscCall(KSPSetTolerances(ksp,1e-8,1e-35,PETSC_DEFAULT,PETSC_DEFAULT));
115: PetscCall(STView(st,NULL));
116: PetscCall(EPSGetDS(eps,&ds));
117: PetscCall(DSView(ds,NULL));
119: PetscCall(EPSSetFromOptions(eps));
120: PetscCall(EPSSolve(eps));
121: PetscCall(EPSGetConvergedReason(eps,&reason));
122: PetscCall(EPSGetIterationNumber(eps,&its));
123: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Finished - converged reason = %d, its=%" PetscInt_FMT "\n",(int)reason,its));
125: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126: Display solution and clean up
127: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
129: PetscCall(EPSDestroy(&eps));
130: PetscCall(MatDestroy(&A));
131: PetscCall(SlepcFinalize());
132: return 0;
133: }
135: /*TEST
137: test:
138: suffix: 1
139: args: -eps_ncv 14
140: filter: sed -e "s/00001/00000/" | sed -e "s/4.99999/5.00000/" | sed -e "s/5.99999/6.00000/"
142: TEST*/