Actual source code: test13.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2012, 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[] = "Test EPSSetArbitrarySelection.\n\n";
24: #include <slepceps.h>
28: PetscErrorCode MyArbitrarySelection(PetscScalar eigr,PetscScalar eigi,Vec xr,Vec xi,PetscScalar *rr,PetscScalar *ri,void *ctx)
29: {
30: PetscErrorCode ierr;
31: Vec xref = *(Vec*)ctx;
33: VecDot(xr,xref,rr);
34: *rr = PetscAbsScalar(*rr);
35: if (ri) *ri = 0.0;
36: return(0);
37: }
41: int main(int argc,char **argv)
42: {
43: Mat A; /* problem matrices */
44: EPS eps; /* eigenproblem solver context */
45: PetscScalar seigr,seigi,value[3];
46: PetscReal tol=1000*PETSC_MACHINE_EPSILON;
47: Vec sxr,sxi;
48: PetscInt n=30,i,Istart,Iend,col[3],nconv;
49: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
52: SlepcInitialize(&argc,&argv,(char*)0,help);
54: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
55: PetscPrintf(PETSC_COMM_WORLD,"\nTridiagonal with zero diagonal, n=%D\n\n",n);
57: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58: Create matrix tridiag([-1 0 -1])
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60: MatCreate(PETSC_COMM_WORLD,&A);
61: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
62: MatSetFromOptions(A);
63: MatSetUp(A);
64:
65: MatGetOwnershipRange(A,&Istart,&Iend);
66: if (Istart==0) FirstBlock=PETSC_TRUE;
67: if (Iend==n) LastBlock=PETSC_TRUE;
68: value[0]=-1.0; value[1]=0.0; value[2]=-1.0;
69: for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) {
70: col[0]=i-1; col[1]=i; col[2]=i+1;
71: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
72: }
73: if (LastBlock) {
74: i=n-1; col[0]=n-2; col[1]=n-1;
75: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
76: }
77: if (FirstBlock) {
78: i=0; col[0]=0; col[1]=1; value[0]=0.0; value[1]=-1.0;
79: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
80: }
82: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
83: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
85: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86: Create the eigensolver
87: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
88: EPSCreate(PETSC_COMM_WORLD,&eps);
89: EPSSetProblemType(eps,EPS_HEP);
90: EPSSetTolerances(eps,tol,PETSC_DECIDE);
91: EPSSetOperators(eps,A,PETSC_NULL);
92: EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
93: EPSSetFromOptions(eps);
95: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96: Solve eigenproblem and store some solution
97: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98: EPSSolve(eps);
99: MatGetVecs(A,&sxr,PETSC_NULL);
100: MatGetVecs(A,&sxi,PETSC_NULL);
101: EPSGetConverged(eps,&nconv);
102: if (nconv>0) {
103: EPSGetEigenpair(eps,0,&seigr,&seigi,sxr,sxi);
104: EPSPrintSolution(eps,PETSC_NULL);
106: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107: Solve eigenproblem using an arbitrary selection
108: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109: EPSSetArbitrarySelection(eps,MyArbitrarySelection,&sxr);
110: EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE);
111: EPSSolve(eps);
112: EPSPrintSolution(eps,PETSC_NULL);
113: } else {
114: PetscPrintf(PETSC_COMM_WORLD,"Problem: no eigenpairs converged.\n");
115: }
116:
117: EPSDestroy(&eps);
118: VecDestroy(&sxr);
119: VecDestroy(&sxi);
120: MatDestroy(&A);
121: SlepcFinalize();
122: return 0;
123: }