Actual source code: test2.c
slepc-3.7.1 2016-05-27
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, 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[] = "Tests multiple calls to EPSSolve with the same matrix.\n\n";
24: #include <slepceps.h>
28: int main(int argc,char **argv)
29: {
30: Mat A; /* problem matrix */
31: EPS eps; /* eigenproblem solver context */
32: ST st;
33: PetscReal tol=PetscMax(1000*PETSC_MACHINE_EPSILON,1e-9);
34: PetscInt n=30,i,Istart,Iend;
35: PetscBool flg;
38: SlepcInitialize(&argc,&argv,(char*)0,help);
40: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
41: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian Eigenproblem, n=%D\n\n",n);
43: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
44: Compute the operator matrix that defines the eigensystem, Ax=kx
45: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
47: MatCreate(PETSC_COMM_WORLD,&A);
48: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
49: MatSetFromOptions(A);
50: MatSetUp(A);
52: MatGetOwnershipRange(A,&Istart,&Iend);
53: for (i=Istart;i<Iend;i++) {
54: if (i>0) { MatSetValue(A,i,i-1,-1.0,INSERT_VALUES); }
55: if (i<n-1) { MatSetValue(A,i,i+1,-1.0,INSERT_VALUES); }
56: MatSetValue(A,i,i,2.0,INSERT_VALUES);
57: }
58: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
59: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
61: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62: Create the eigensolver
63: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
64: EPSCreate(PETSC_COMM_WORLD,&eps);
65: EPSSetOperators(eps,A,NULL);
66: EPSSetProblemType(eps,EPS_HEP);
67: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
68: EPSSetFromOptions(eps);
70: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71: Solve for largest eigenvalues
72: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73: EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
74: EPSSolve(eps);
75: PetscPrintf(PETSC_COMM_WORLD," - - - Largest eigenvalues - - -\n");
76: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
78: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79: Solve for smallest eigenvalues
80: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
81: EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
82: EPSSolve(eps);
83: PetscPrintf(PETSC_COMM_WORLD," - - - Smallest eigenvalues - - -\n");
84: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
86: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87: Solve for interior eigenvalues (target=2.1)
88: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
89: EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);
90: EPSSetTarget(eps,2.1);
91: PetscObjectTypeCompare((PetscObject)eps,EPSLANCZOS,&flg);
92: if (flg) {
93: EPSGetST(eps,&st);
94: STSetType(st,STSINVERT);
95: } else {
96: PetscObjectTypeCompare((PetscObject)eps,EPSKRYLOVSCHUR,&flg);
97: if (!flg) {
98: PetscObjectTypeCompare((PetscObject)eps,EPSARNOLDI,&flg);
99: }
100: if (flg) {
101: EPSSetExtraction(eps,EPS_HARMONIC);
102: }
103: }
104: EPSSolve(eps);
105: PetscPrintf(PETSC_COMM_WORLD," - - - Interior eigenvalues - - -\n");
106: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
108: EPSDestroy(&eps);
109: MatDestroy(&A);
110: SlepcFinalize();
111: return ierr;
112: }