Actual source code: test14.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[] = "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,NULL,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(MatGetOwnershipRange(A,&Istart,&Iend));
44: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(A,i,i,i+1,INSERT_VALUES));
45: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
46: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
48: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49: Create eigensolver and test interface functions
50: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
51: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
52: PetscCall(EPSSetOperators(eps,A,NULL));
53: PetscCall(EPSGetOperators(eps,&B,NULL));
54: PetscCall(MatView(B,NULL));
56: PetscCall(EPSSetType(eps,EPSKRYLOVSCHUR));
57: PetscCall(EPSGetType(eps,&type));
58: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Type set to %s\n",type));
60: PetscCall(EPSGetProblemType(eps,&ptype));
61: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Problem type before changing = %d",(int)ptype));
62: PetscCall(EPSSetProblemType(eps,EPS_HEP));
63: PetscCall(EPSGetProblemType(eps,&ptype));
64: PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %d.",(int)ptype));
65: PetscCall(EPSIsGeneralized(eps,&flg));
66: if (flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD," generalized"));
67: PetscCall(EPSIsHermitian(eps,&flg));
68: if (flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD," hermitian"));
69: PetscCall(EPSIsPositive(eps,&flg));
70: if (flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD," positive"));
72: PetscCall(EPSGetExtraction(eps,&extr));
73: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n Extraction before changing = %d",(int)extr));
74: PetscCall(EPSSetExtraction(eps,EPS_HARMONIC));
75: PetscCall(EPSGetExtraction(eps,&extr));
76: PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)extr));
78: PetscCall(EPSSetBalance(eps,EPS_BALANCE_ONESIDE,8,1e-6));
79: PetscCall(EPSGetBalance(eps,&bal,&its,&cut));
80: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Balance: %s, its=%" PetscInt_FMT ", cutoff=%g\n",EPSBalanceTypes[bal],its,(double)cut));
82: PetscCall(EPSSetPurify(eps,PETSC_FALSE));
83: PetscCall(EPSGetPurify(eps,&pur));
84: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Eigenvector purification: %s\n",pur?"on":"off"));
85: PetscCall(EPSGetTrackAll(eps,&track));
87: PetscCall(EPSSetTarget(eps,4.8));
88: PetscCall(EPSGetTarget(eps,&target));
89: PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
90: PetscCall(EPSGetWhichEigenpairs(eps,&which));
91: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Which = %d, target = %g\n",(int)which,(double)PetscRealPart(target)));
93: PetscCall(EPSSetDimensions(eps,4,PETSC_DETERMINE,PETSC_DETERMINE));
94: PetscCall(EPSGetDimensions(eps,&nev,&ncv,&mpd));
95: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Dimensions: nev=%" PetscInt_FMT ", ncv=%" PetscInt_FMT ", mpd=%" PetscInt_FMT "\n",nev,ncv,mpd));
97: PetscCall(EPSSetTolerances(eps,2.2e-4,200));
98: PetscCall(EPSGetTolerances(eps,&tol,&its));
99: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Tolerance = %.5f, max_its = %" PetscInt_FMT "\n",(double)tol,its));
101: PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_ABS));
102: PetscCall(EPSGetConvergenceTest(eps,&conv));
103: PetscCall(EPSSetStoppingTest(eps,EPS_STOP_BASIC));
104: PetscCall(EPSGetStoppingTest(eps,&stop));
105: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Convergence test = %d, stopping test = %d\n",(int)conv,(int)stop));
107: PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf));
108: PetscCall(EPSMonitorSet(eps,(PetscErrorCode (*)(EPS,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))EPSMonitorFirst,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy));
109: PetscCall(EPSMonitorCancel(eps));
111: PetscCall(EPSGetST(eps,&st));
112: PetscCall(STGetKSP(st,&ksp));
113: PetscCall(KSPSetTolerances(ksp,1e-8,1e-35,PETSC_CURRENT,PETSC_CURRENT));
114: PetscCall(STView(st,NULL));
115: PetscCall(EPSGetDS(eps,&ds));
116: PetscCall(DSView(ds,NULL));
118: PetscCall(EPSSetFromOptions(eps));
119: PetscCall(EPSSolve(eps));
120: PetscCall(EPSGetConvergedReason(eps,&reason));
121: PetscCall(EPSGetIterationNumber(eps,&its));
122: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Finished - converged reason = %d, its=%" PetscInt_FMT "\n",(int)reason,its));
124: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125: Display solution and clean up
126: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
128: PetscCall(EPSDestroy(&eps));
129: PetscCall(MatDestroy(&A));
130: PetscCall(SlepcFinalize());
131: return 0;
132: }
134: /*TEST
136: test:
137: suffix: 1
138: args: -eps_ncv 14
139: filter: sed -e "s/00001/00000/" | sed -e "s/4.99999/5.00000/" | sed -e "s/5.99999/6.00000/"
141: TEST*/