Actual source code: test2.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7:
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=1000*PETSC_MACHINE_EPSILON;
34: PetscScalar value[3];
35: PetscInt n=30,i,Istart,Iend,col[3];
36: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE,flg;
39: SlepcInitialize(&argc,&argv,(char*)0,help);
41: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
42: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian Eigenproblem, n=%D\n\n",n);
44: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
45: Compute the operator matrix that defines the eigensystem, Ax=kx
46: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
48: MatCreate(PETSC_COMM_WORLD,&A);
49: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
50: MatSetFromOptions(A);
51:
52: MatGetOwnershipRange(A,&Istart,&Iend);
53: if (Istart==0) FirstBlock=PETSC_TRUE;
54: if (Iend==n) LastBlock=PETSC_TRUE;
55: value[0]=-1.0; value[1]=2.0; value[2]=-1.0;
56: for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) {
57: col[0]=i-1; col[1]=i; col[2]=i+1;
58: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
59: }
60: if (LastBlock) {
61: i=n-1; col[0]=n-2; col[1]=n-1;
62: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
63: }
64: if (FirstBlock) {
65: i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0;
66: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
67: }
69: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
70: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
72: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73: Create the eigensolver
74: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
75: EPSCreate(PETSC_COMM_WORLD,&eps);
76: EPSSetOperators(eps,A,PETSC_NULL);
77: EPSSetProblemType(eps,EPS_HEP);
78: EPSSetTolerances(eps,tol,PETSC_DECIDE);
79: EPSSetFromOptions(eps);
81: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82: Solve for largest eigenvalues
83: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
84: EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
85: EPSSolve(eps);
86: PetscPrintf(PETSC_COMM_WORLD," - - - Largest eigenvalues - - -\n");
87: EPSPrintSolution(eps,PETSC_NULL);
89: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
90: Solve for smallest eigenvalues
91: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
92: EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
93: EPSSolve(eps);
94: PetscPrintf(PETSC_COMM_WORLD," - - - Smallest eigenvalues - - -\n");
95: EPSPrintSolution(eps,PETSC_NULL);
96:
97: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98: Solve for interior eigenvalues (target=2.1)
99: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100: EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);
101: EPSSetTarget(eps,2.1);
102: PetscTypeCompare((PetscObject)eps,EPSLANCZOS,&flg);
103: if (flg) {
104: EPSGetST(eps,&st);
105: STSetType(st,STSINVERT);
106: } else {
107: PetscTypeCompare((PetscObject)eps,EPSKRYLOVSCHUR,&flg);
108: if (!flg) {
109: PetscTypeCompare((PetscObject)eps,EPSARNOLDI,&flg);
110: }
111: if (flg) {
112: EPSSetExtraction(eps,EPS_HARMONIC);
113: }
114: }
115: EPSSolve(eps);
116: PetscPrintf(PETSC_COMM_WORLD," - - - Interior eigenvalues - - -\n");
117: EPSPrintSolution(eps,PETSC_NULL);
118:
119: EPSDestroy(&eps);
120: MatDestroy(&A);
121: SlepcFinalize();
122: return 0;
123: }