Actual source code: test14.c
slepc-3.5.0 2014-07-29
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2014, 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: DS ds;
34: PetscReal cut,tol;
35: PetscScalar target;
36: PetscInt n=20,i,its,nev,ncv,mpd,Istart,Iend;
37: PetscBool flg;
38: EPSConvergedReason reason;
39: EPSType type;
40: EPSExtraction extr;
41: EPSBalance bal;
42: EPSWhich which;
43: EPSConv conv;
44: EPSProblemType ptype;
45: PetscErrorCode ierr;
47: SlepcInitialize(&argc,&argv,(char*)0,help);
48: PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Eigenproblem, n=%D\n\n",n);
50: MatCreate(PETSC_COMM_WORLD,&A);
51: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
52: MatSetFromOptions(A);
53: MatSetUp(A);
54: MatGetOwnershipRange(A,&Istart,&Iend);
55: for (i=Istart;i<Iend;i++) {
56: MatSetValue(A,i,i,i+1,INSERT_VALUES);
57: }
58: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
59: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
61: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62: Create eigensolver and test interface functions
63: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
64: EPSCreate(PETSC_COMM_WORLD,&eps);
65: EPSSetOperators(eps,A,NULL);
66: EPSGetOperators(eps,&B,NULL);
67: MatView(B,NULL);
69: EPSSetType(eps,EPSKRYLOVSCHUR);
70: EPSGetType(eps,&type);
71: PetscPrintf(PETSC_COMM_WORLD," Type set to %s\n",type);
73: EPSGetProblemType(eps,&ptype);
74: PetscPrintf(PETSC_COMM_WORLD," Problem type before changing = %D",ptype);
75: EPSSetProblemType(eps,EPS_HEP);
76: EPSGetProblemType(eps,&ptype);
77: PetscPrintf(PETSC_COMM_WORLD," ... changed to %D.",ptype);
78: EPSIsGeneralized(eps,&flg);
79: if (flg) { PetscPrintf(PETSC_COMM_WORLD," generalized"); }
80: EPSIsHermitian(eps,&flg);
81: if (flg) { PetscPrintf(PETSC_COMM_WORLD," hermitian"); }
82: EPSIsPositive(eps,&flg);
83: if (flg) { PetscPrintf(PETSC_COMM_WORLD," positive"); }
85: EPSGetExtraction(eps,&extr);
86: PetscPrintf(PETSC_COMM_WORLD,"\n Extraction before changing = %D",extr);
87: EPSSetExtraction(eps,EPS_HARMONIC);
88: EPSGetExtraction(eps,&extr);
89: PetscPrintf(PETSC_COMM_WORLD," ... changed to %D\n",extr);
91: EPSSetBalance(eps,EPS_BALANCE_ONESIDE,8,1e-6);
92: EPSGetBalance(eps,&bal,&its,&cut);
93: PetscPrintf(PETSC_COMM_WORLD," Balance: %s, its=%D, cutoff=%g\n",EPSBalanceTypes[bal],its,(double)cut);
95: EPSSetTarget(eps,4.8);
96: EPSGetTarget(eps,&target);
97: EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);
98: EPSGetWhichEigenpairs(eps,&which);
99: PetscPrintf(PETSC_COMM_WORLD," Which = %D, target = %g\n",which,(double)PetscRealPart(target));
101: EPSSetDimensions(eps,4,PETSC_DEFAULT,PETSC_DEFAULT);
102: EPSGetDimensions(eps,&nev,&ncv,&mpd);
103: PetscPrintf(PETSC_COMM_WORLD," Dimensions: nev=%D, ncv=%D, mpd=%D\n",nev,ncv,mpd);
105: EPSSetTolerances(eps,2.2e-4,200);
106: EPSGetTolerances(eps,&tol,&its);
107: PetscPrintf(PETSC_COMM_WORLD," Tolerance = %.5f, max_its = %D\n",(double)tol,its);
109: EPSSetConvergenceTest(eps,EPS_CONV_ABS);
110: EPSGetConvergenceTest(eps,&conv);
111: PetscPrintf(PETSC_COMM_WORLD," Convergence test = %D\n",conv);
113: EPSMonitorSet(eps,EPSMonitorFirst,NULL,NULL);
114: EPSMonitorCancel(eps);
116: EPSGetST(eps,&st);
117: STView(st,NULL);
118: EPSGetDS(eps,&ds);
119: DSView(ds,NULL);
121: EPSSetFromOptions(eps);
122: EPSSolve(eps);
123: EPSGetConvergedReason(eps,&reason);
124: EPSGetIterationNumber(eps,&its);
125: PetscPrintf(PETSC_COMM_WORLD," Finished - converged reason = %D, its=%D\n",reason,its);
127: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128: Display solution and clean up
129: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130: EPSPrintSolution(eps,NULL);
131: EPSDestroy(&eps);
132: MatDestroy(&A);
133: SlepcFinalize();
134: return 0;
135: }