Actual source code: test12.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[] = "Diagonal eigenproblem.\n\n"
23: "The command line options are:\n"
24: " -n <n>, where <n> = number of grid subdivisions = matrix dimension.\n"
25: " -seed <s>, where <s> = seed for random number generation.\n\n";
27: #include <slepceps.h>
31: PetscErrorCode MyPCApply(PC pc,Vec x,Vec y)
32: {
35: VecCopy(x,y);
36: return(0);
37: }
41: int main(int argc,char **argv)
42: {
43: Mat A; /* problem matrix */
44: EPS eps; /* eigenproblem solver context */
45: Vec v0; /* initial vector */
46: PetscRandom rand;
47: PetscReal tol=1000*PETSC_MACHINE_EPSILON;
48: PetscInt n=30,i,Istart,Iend,seed=0x12345678;
50: ST st;
51: KSP ksp;
52: PC pc;
54: SlepcInitialize(&argc,&argv,(char*)0,help);
56: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
57: PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Eigenproblem, n=%d\n\n",n);
59: MatCreate(PETSC_COMM_WORLD,&A);
60: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
61: MatSetFromOptions(A);
62: MatSetUp(A);
63: MatGetOwnershipRange(A,&Istart,&Iend);
64: for (i=Istart;i<Iend;i++) {
65: MatSetValue(A,i,i,i+1,INSERT_VALUES);
66: }
67: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
68: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
70: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71: Solve the eigensystem
72: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73: EPSCreate(PETSC_COMM_WORLD,&eps);
74: EPSSetOperators(eps,A,PETSC_NULL);
75: EPSSetProblemType(eps,EPS_HEP);
76: EPSSetTolerances(eps,tol,PETSC_DECIDE);
77: EPSSetFromOptions(eps);
78: EPSGetST(eps,&st);
79: STGetKSP(st,&ksp);
80: KSPGetPC(ksp,&pc);
81: PCSetType(pc,PCSHELL);
82: PCShellSetApply(pc,MyPCApply);
83: STPrecondSetMatForPC(st,A);
85: /* set random initial vector */
86: MatGetVecs(A,&v0,PETSC_NULL);
87: PetscRandomCreate(PETSC_COMM_WORLD,&rand);
88: PetscRandomSetFromOptions(rand);
89: PetscOptionsGetInt(PETSC_NULL,"-seed",&seed,PETSC_NULL);
90: PetscRandomSetSeed(rand,seed);
91: PetscRandomSeed(rand);
92: VecSetRandom(v0,rand);
93: EPSSetInitialSpace(eps,1,&v0);
94: /* call the solver */
95: EPSSolve(eps);
97: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98: Display solution and clean up
99: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100: EPSPrintSolution(eps,PETSC_NULL);
101: EPSDestroy(&eps);
102: MatDestroy(&A);
103: VecDestroy(&v0);
104: PetscRandomDestroy(&rand);
105: SlepcFinalize();
106: return 0;
107: }