Actual source code: test6.c
slepc-3.7.2 2016-07-19
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[] = "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: int main(int argc,char **argv)
32: {
33: Mat A; /* problem matrix */
34: EPS eps; /* eigenproblem solver context */
35: Vec v0; /* initial vector */
36: PetscRandom rand;
37: PetscReal tol=1000*PETSC_MACHINE_EPSILON;
38: PetscInt n=30,i,Istart,Iend,seed=0x12345678;
41: SlepcInitialize(&argc,&argv,(char*)0,help);
43: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
44: PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Eigenproblem, n=%D\n\n",n);
46: MatCreate(PETSC_COMM_WORLD,&A);
47: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
48: MatSetFromOptions(A);
49: MatSetUp(A);
50: MatGetOwnershipRange(A,&Istart,&Iend);
51: for (i=Istart;i<Iend;i++) {
52: MatSetValue(A,i,i,i+1,INSERT_VALUES);
53: }
54: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
55: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
57: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58: Solve the eigensystem
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60: EPSCreate(PETSC_COMM_WORLD,&eps);
61: EPSSetOperators(eps,A,NULL);
62: EPSSetProblemType(eps,EPS_HEP);
63: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
64: EPSSetFromOptions(eps);
65: /* set random initial vector */
66: MatCreateVecs(A,&v0,NULL);
67: PetscRandomCreate(PETSC_COMM_WORLD,&rand);
68: PetscRandomSetFromOptions(rand);
69: PetscOptionsGetInt(NULL,NULL,"-seed",&seed,NULL);
70: PetscRandomSetSeed(rand,seed);
71: PetscRandomSeed(rand);
72: VecSetRandom(v0,rand);
73: EPSSetInitialSpace(eps,1,&v0);
74: /* call the solver */
75: EPSSolve(eps);
77: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
78: Display solution and clean up
79: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
80: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
81: EPSDestroy(&eps);
82: MatDestroy(&A);
83: VecDestroy(&v0);
84: PetscRandomDestroy(&rand);
85: SlepcFinalize();
86: return ierr;
87: }