Actual source code: ex2.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2011, 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[] = "Standard symmetric eigenproblem corresponding to the Laplacian operator in 2 dimensions.\n\n"
23: "The command line options are:\n"
24: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
25: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";
27: #include <slepceps.h>
31: int main(int argc,char **argv)
32: {
33: Mat A; /* operator matrix */
34: EPS eps; /* eigenproblem solver context */
35: const EPSType type;
36: PetscReal tol;
37: PetscInt N,n=10,m,Istart,Iend,II,nev,maxit,i,j,its;
38: PetscBool flag;
41: SlepcInitialize(&argc,&argv,(char*)0,help);
43: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
44: PetscOptionsGetInt(PETSC_NULL,"-m",&m,&flag);
45: if(!flag) m=n;
46: N = n*m;
47: PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem, N=%D (%Dx%D grid)\n\n",N,n,m);
49: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50: Compute the operator matrix that defines the eigensystem, Ax=kx
51: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
53: MatCreate(PETSC_COMM_WORLD,&A);
54: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
55: MatSetFromOptions(A);
56:
57: MatGetOwnershipRange(A,&Istart,&Iend);
58: for (II=Istart;II<Iend;II++) {
59: i = II/n; j = II-i*n;
60: if(i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
61: if(i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
62: if(j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
63: if(j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
64: MatSetValue(A,II,II,4.0,INSERT_VALUES);
65: }
67: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
68: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
70: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71: Create the eigensolver and set various options
72: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74: /*
75: Create eigensolver context
76: */
77: EPSCreate(PETSC_COMM_WORLD,&eps);
79: /*
80: Set operators. In this case, it is a standard eigenvalue problem
81: */
82: EPSSetOperators(eps,A,PETSC_NULL);
83: EPSSetProblemType(eps,EPS_HEP);
85: /*
86: Set solver parameters at runtime
87: */
88: EPSSetFromOptions(eps);
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Solve the eigensystem
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94: EPSSolve(eps);
95: EPSGetIterationNumber(eps,&its);
96: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",its);
98: /*
99: Optional: Get some information from the solver and display it
100: */
101: EPSGetType(eps,&type);
102: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
103: EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);
104: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
105: EPSGetTolerances(eps,&tol,&maxit);
106: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4G, maxit=%D\n",tol,maxit);
108: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109: Display solution and clean up
110: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
112: EPSPrintSolution(eps,PETSC_NULL);
113: EPSDestroy(&eps);
114: MatDestroy(&A);
115: SlepcFinalize();
116: return 0;
117: }