Actual source code: test2.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2013, 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=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(NULL,"-n",&n,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: MatSetUp(A);
53: MatGetOwnershipRange(A,&Istart,&Iend);
54: if (Istart==0) FirstBlock=PETSC_TRUE;
55: if (Iend==n) LastBlock=PETSC_TRUE;
56: value[0]=-1.0; value[1]=2.0; value[2]=-1.0;
57: for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) {
58: col[0]=i-1; col[1]=i; col[2]=i+1;
59: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
60: }
61: if (LastBlock) {
62: i=n-1; col[0]=n-2; col[1]=n-1;
63: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
64: }
65: if (FirstBlock) {
66: i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0;
67: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
68: }
70: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
71: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
73: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74: Create the eigensolver
75: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
76: EPSCreate(PETSC_COMM_WORLD,&eps);
77: EPSSetOperators(eps,A,NULL);
78: EPSSetProblemType(eps,EPS_HEP);
79: EPSSetTolerances(eps,tol,PETSC_DECIDE);
80: EPSSetFromOptions(eps);
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Solve for largest eigenvalues
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85: EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
86: EPSSolve(eps);
87: PetscPrintf(PETSC_COMM_WORLD," - - - Largest eigenvalues - - -\n");
88: EPSPrintSolution(eps,NULL);
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Solve for smallest eigenvalues
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93: EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
94: EPSSolve(eps);
95: PetscPrintf(PETSC_COMM_WORLD," - - - Smallest eigenvalues - - -\n");
96: EPSPrintSolution(eps,NULL);
98: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99: Solve for interior eigenvalues (target=2.1)
100: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
101: EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);
102: EPSSetTarget(eps,2.1);
103: PetscObjectTypeCompare((PetscObject)eps,EPSLANCZOS,&flg);
104: if (flg) {
105: EPSGetST(eps,&st);
106: STSetType(st,STSINVERT);
107: } else {
108: PetscObjectTypeCompare((PetscObject)eps,EPSKRYLOVSCHUR,&flg);
109: if (!flg) {
110: PetscObjectTypeCompare((PetscObject)eps,EPSARNOLDI,&flg);
111: }
112: if (flg) {
113: EPSSetExtraction(eps,EPS_HARMONIC);
114: }
115: }
116: EPSSolve(eps);
117: PetscPrintf(PETSC_COMM_WORLD," - - - Interior eigenvalues - - -\n");
118: EPSPrintSolution(eps,NULL);
120: EPSDestroy(&eps);
121: MatDestroy(&A);
122: SlepcFinalize();
123: return 0;
124: }