Actual source code: ex1.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2013, 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[] = "Standard symmetric eigenproblem corresponding to the Laplacian operator in 1 dimension.\n\n"
23: "The command line options are:\n"
24: " -n <n>, where <n> = number of grid subdivisions = matrix dimension.\n\n";
26: #include <slepceps.h>
30: int main(int argc,char **argv)
31: {
32: Mat A; /* problem matrix */
33: EPS eps; /* eigenproblem solver context */
34: EPSType type;
35: PetscReal error,tol,re,im;
36: PetscScalar kr,ki,value[3];
37: Vec xr,xi;
38: PetscInt n=30,i,Istart,Iend,col[3],nev,maxit,its,nconv;
39: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
42: SlepcInitialize(&argc,&argv,(char*)0,help);
44: PetscOptionsGetInt(NULL,"-n",&n,NULL);
45: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian Eigenproblem, n=%D\n\n",n);
47: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48: Compute the operator matrix that defines the eigensystem, Ax=kx
49: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
51: MatCreate(PETSC_COMM_WORLD,&A);
52: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
53: MatSetFromOptions(A);
54: MatSetUp(A);
56: MatGetOwnershipRange(A,&Istart,&Iend);
57: if (Istart==0) FirstBlock=PETSC_TRUE;
58: if (Iend==n) LastBlock=PETSC_TRUE;
59: value[0]=-1.0; value[1]=2.0; value[2]=-1.0;
60: for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) {
61: col[0]=i-1; col[1]=i; col[2]=i+1;
62: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
63: }
64: if (LastBlock) {
65: i=n-1; col[0]=n-2; col[1]=n-1;
66: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
67: }
68: if (FirstBlock) {
69: i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0;
70: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
71: }
73: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
74: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
76: MatGetVecs(A,NULL,&xr);
77: MatGetVecs(A,NULL,&xi);
79: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80: Create the eigensolver and set various options
81: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82: /*
83: Create eigensolver context
84: */
85: EPSCreate(PETSC_COMM_WORLD,&eps);
87: /*
88: Set operators. In this case, it is a standard eigenvalue problem
89: */
90: EPSSetOperators(eps,A,NULL);
91: EPSSetProblemType(eps,EPS_HEP);
93: /*
94: Set solver parameters at runtime
95: */
96: EPSSetFromOptions(eps);
98: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99: Solve the eigensystem
100: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
102: EPSSolve(eps);
103: /*
104: Optional: Get some information from the solver and display it
105: */
106: EPSGetIterationNumber(eps,&its);
107: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",its);
108: EPSGetType(eps,&type);
109: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
110: EPSGetDimensions(eps,&nev,NULL,NULL);
111: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
112: EPSGetTolerances(eps,&tol,&maxit);
113: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4G, maxit=%D\n",tol,maxit);
115: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116: Display solution and clean up
117: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118: /*
119: Get number of converged approximate eigenpairs
120: */
121: EPSGetConverged(eps,&nconv);
122: PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %D\n\n",nconv);
124: if (nconv>0) {
125: /*
126: Display eigenvalues and relative errors
127: */
128: PetscPrintf(PETSC_COMM_WORLD,
129: " k ||Ax-kx||/||kx||\n"
130: " ----------------- ------------------\n");
132: for (i=0;i<nconv;i++) {
133: /*
134: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
135: ki (imaginary part)
136: */
137: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
138: /*
139: Compute the relative error associated to each eigenpair
140: */
141: EPSComputeRelativeError(eps,i,&error);
143: #if defined(PETSC_USE_COMPLEX)
144: re = PetscRealPart(kr);
145: im = PetscImaginaryPart(kr);
146: #else
147: re = kr;
148: im = ki;
149: #endif
150: if (im!=0.0) {
151: PetscPrintf(PETSC_COMM_WORLD," %9F%+9F j %12G\n",re,im,error);
152: } else {
153: PetscPrintf(PETSC_COMM_WORLD," %12F %12G\n",re,error);
154: }
155: }
156: PetscPrintf(PETSC_COMM_WORLD,"\n");
157: }
159: /*
160: Free work space
161: */
162: EPSDestroy(&eps);
163: MatDestroy(&A);
164: VecDestroy(&xr);
165: VecDestroy(&xi);
166: SlepcFinalize();
167: return 0;
168: }