Actual source code: test13.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[] = "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;
34: VecDot(xr,xref,rr);
35: *rr = PetscAbsScalar(*rr);
36: if (ri) *ri = 0.0;
37: return(0);
38: }
42: int main(int argc,char **argv)
43: {
44: Mat A; /* problem matrices */
45: EPS eps; /* eigenproblem solver context */
46: PetscScalar seigr,seigi;
47: PetscReal tol=1000*PETSC_MACHINE_EPSILON;
48: Vec sxr,sxi;
49: PetscInt n=30,i,Istart,Iend,nconv;
52: SlepcInitialize(&argc,&argv,(char*)0,help);
54: PetscOptionsGetInt(NULL,NULL,"-n",&n,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);
65: MatGetOwnershipRange(A,&Istart,&Iend);
66: for (i=Istart;i<Iend;i++) {
67: if (i>0) { MatSetValue(A,i,i-1,-1.0,INSERT_VALUES); }
68: if (i<n-1) { MatSetValue(A,i,i+1,-1.0,INSERT_VALUES); }
69: }
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: EPSSetProblemType(eps,EPS_HEP);
78: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
79: EPSSetOperators(eps,A,NULL);
80: EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
81: EPSSetFromOptions(eps);
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Solve eigenproblem and store some solution
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: EPSSolve(eps);
87: MatCreateVecs(A,&sxr,NULL);
88: MatCreateVecs(A,&sxi,NULL);
89: EPSGetConverged(eps,&nconv);
90: if (nconv>0) {
91: EPSGetEigenpair(eps,0,&seigr,&seigi,sxr,sxi);
92: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
94: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: Solve eigenproblem using an arbitrary selection
96: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97: EPSSetArbitrarySelection(eps,MyArbitrarySelection,&sxr);
98: EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE);
99: EPSSolve(eps);
100: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
101: } else {
102: PetscPrintf(PETSC_COMM_WORLD,"Problem: no eigenpairs converged.\n");
103: }
105: EPSDestroy(&eps);
106: VecDestroy(&sxr);
107: VecDestroy(&sxi);
108: MatDestroy(&A);
109: SlepcFinalize();
110: return ierr;
111: }