Actual source code: test14.c
slepc-3.6.1 2015-09-03
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: static char help[] = "Test EPS interface functions.\n\n";
24: #include <slepceps.h>
28: int main(int argc,char **argv)
29: {
30: Mat A,B; /* problem matrix */
31: EPS eps; /* eigenproblem solver context */
32: ST st;
33: KSP ksp;
34: DS ds;
35: PetscReal cut,tol;
36: PetscScalar target;
37: PetscInt n=20,i,its,nev,ncv,mpd,Istart,Iend;
38: PetscBool flg;
39: EPSConvergedReason reason;
40: EPSType type;
41: EPSExtraction extr;
42: EPSBalance bal;
43: EPSWhich which;
44: EPSConv conv;
45: EPSProblemType ptype;
46: PetscErrorCode ierr;
48: SlepcInitialize(&argc,&argv,(char*)0,help);
49: PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Eigenproblem, n=%D\n\n",n);
51: MatCreate(PETSC_COMM_WORLD,&A);
52: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
53: MatSetFromOptions(A);
54: MatSetUp(A);
55: MatGetOwnershipRange(A,&Istart,&Iend);
56: for (i=Istart;i<Iend;i++) {
57: MatSetValue(A,i,i,i+1,INSERT_VALUES);
58: }
59: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
60: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
62: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63: Create eigensolver and test interface functions
64: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65: EPSCreate(PETSC_COMM_WORLD,&eps);
66: EPSSetOperators(eps,A,NULL);
67: EPSGetOperators(eps,&B,NULL);
68: MatView(B,NULL);
70: EPSSetType(eps,EPSKRYLOVSCHUR);
71: EPSGetType(eps,&type);
72: PetscPrintf(PETSC_COMM_WORLD," Type set to %s\n",type);
74: EPSGetProblemType(eps,&ptype);
75: PetscPrintf(PETSC_COMM_WORLD," Problem type before changing = %d",(int)ptype);
76: EPSSetProblemType(eps,EPS_HEP);
77: EPSGetProblemType(eps,&ptype);
78: PetscPrintf(PETSC_COMM_WORLD," ... changed to %d.",(int)ptype);
79: EPSIsGeneralized(eps,&flg);
80: if (flg) { PetscPrintf(PETSC_COMM_WORLD," generalized"); }
81: EPSIsHermitian(eps,&flg);
82: if (flg) { PetscPrintf(PETSC_COMM_WORLD," hermitian"); }
83: EPSIsPositive(eps,&flg);
84: if (flg) { PetscPrintf(PETSC_COMM_WORLD," positive"); }
86: EPSGetExtraction(eps,&extr);
87: PetscPrintf(PETSC_COMM_WORLD,"\n Extraction before changing = %d",(int)extr);
88: EPSSetExtraction(eps,EPS_HARMONIC);
89: EPSGetExtraction(eps,&extr);
90: PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)extr);
92: EPSSetBalance(eps,EPS_BALANCE_ONESIDE,8,1e-6);
93: EPSGetBalance(eps,&bal,&its,&cut);
94: PetscPrintf(PETSC_COMM_WORLD," Balance: %s, its=%D, cutoff=%g\n",EPSBalanceTypes[bal],its,(double)cut);
96: EPSSetTarget(eps,4.8);
97: EPSGetTarget(eps,&target);
98: EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);
99: EPSGetWhichEigenpairs(eps,&which);
100: PetscPrintf(PETSC_COMM_WORLD," Which = %d, target = %g\n",(int)which,(double)PetscRealPart(target));
102: EPSSetDimensions(eps,4,PETSC_DEFAULT,PETSC_DEFAULT);
103: EPSGetDimensions(eps,&nev,&ncv,&mpd);
104: PetscPrintf(PETSC_COMM_WORLD," Dimensions: nev=%D, ncv=%D, mpd=%D\n",nev,ncv,mpd);
106: EPSSetTolerances(eps,2.2e-4,200);
107: EPSGetTolerances(eps,&tol,&its);
108: PetscPrintf(PETSC_COMM_WORLD," Tolerance = %.5f, max_its = %D\n",(double)tol,its);
110: EPSSetConvergenceTest(eps,EPS_CONV_ABS);
111: EPSGetConvergenceTest(eps,&conv);
112: PetscPrintf(PETSC_COMM_WORLD," Convergence test = %d\n",(int)conv);
114: EPSMonitorSet(eps,EPSMonitorFirst,NULL,NULL);
115: EPSMonitorCancel(eps);
117: EPSGetST(eps,&st);
118: STGetKSP(st,&ksp);
119: KSPSetTolerances(ksp,1e-8,1e-50,PETSC_DEFAULT,PETSC_DEFAULT);
120: STView(st,NULL);
121: EPSGetDS(eps,&ds);
122: DSView(ds,NULL);
124: EPSSetFromOptions(eps);
125: EPSSolve(eps);
126: EPSGetConvergedReason(eps,&reason);
127: EPSGetIterationNumber(eps,&its);
128: PetscPrintf(PETSC_COMM_WORLD," Finished - converged reason = %d, its=%D\n",(int)reason,its);
130: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131: Display solution and clean up
132: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
134: EPSDestroy(&eps);
135: MatDestroy(&A);
136: SlepcFinalize();
137: return 0;
138: }