Actual source code: ex3.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
6: This file is part of SLEPc. See the README file for conditions of use
7: and additional information.
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: static char help[] = "Solves the same eigenproblem as in example ex2, but using a shell matrix. "
12: "The problem is a standard symmetric eigenproblem corresponding to the 2-D Laplacian operator.\n\n"
13: "The command line options are:\n"
14: " -n <n>, where <n> = number of grid subdivisions in both x and y dimensions.\n\n";
16: #include slepceps.h
17: #include "petscblaslapack.h"
19: /*
20: User-defined routines
21: */
22: PetscErrorCode MatLaplacian2D_Mult( Mat A, Vec x, Vec y );
26: int main( int argc, char **argv )
27: {
28: Mat A; /* operator matrix */
29: EPS eps; /* eigenproblem solver context */
30: EPSType type;
31: PetscReal error, tol, re, im;
32: PetscScalar kr, ki;
33: PetscMPIInt size;
35: PetscInt N, n=10;
36: int nev, maxit, i, its, nconv;
38: SlepcInitialize(&argc,&argv,(char*)0,help);
39: MPI_Comm_size(PETSC_COMM_WORLD,&size);
40: if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
42: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
43: N = n*n;
44: PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem (matrix-free version), N=%d (%dx%d grid)\n\n",N,n,n);
46: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47: Compute the operator matrix that defines the eigensystem, Ax=kx
48: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50: MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,&n,&A);
51: MatSetFromOptions(A);
52: MatShellSetOperation(A,MATOP_MULT,(void(*)())MatLaplacian2D_Mult);
53: MatShellSetOperation(A,MATOP_MULT_TRANSPOSE,(void(*)())MatLaplacian2D_Mult);
55: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
56: Create the eigensolver and set various options
57: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59: /*
60: Create eigensolver context
61: */
62: EPSCreate(PETSC_COMM_WORLD,&eps);
64: /*
65: Set operators. In this case, it is a standard eigenvalue problem
66: */
67: EPSSetOperators(eps,A,PETSC_NULL);
68: EPSSetProblemType(eps,EPS_HEP);
70: /*
71: Set solver parameters at runtime
72: */
73: EPSSetFromOptions(eps);
75: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76: Solve the eigensystem
77: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
79: EPSSolve(eps);
80: EPSGetIterationNumber(eps, &its);
81: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
83: /*
84: Optional: Get some information from the solver and display it
85: */
86: EPSGetType(eps,&type);
87: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
88: EPSGetDimensions(eps,&nev,PETSC_NULL);
89: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
90: EPSGetTolerances(eps,&tol,&maxit);
91: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
93: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94: Display solution and clean up
95: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97: /*
98: Get number of converged approximate eigenpairs
99: */
100: EPSGetConverged(eps,&nconv);
101: PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %d\n\n",nconv);
102:
104: if (nconv>0) {
105: /*
106: Display eigenvalues and relative errors
107: */
108: PetscPrintf(PETSC_COMM_WORLD,
109: " k ||Ax-kx||/||kx||\n"
110: " ----------------- ------------------\n" );
112: for( i=0; i<nconv; i++ ) {
113: /*
114: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
115: ki (imaginary part)
116: */
117: EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
118: /*
119: Compute the relative error associated to each eigenpair
120: */
121: EPSComputeRelativeError(eps,i,&error);
123: #ifdef PETSC_USE_COMPLEX
124: re = PetscRealPart(kr);
125: im = PetscImaginaryPart(kr);
126: #else
127: re = kr;
128: im = ki;
129: #endif
130: if (im!=0.0) {
131: PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12g\n",re,im,error);
132: } else {
133: PetscPrintf(PETSC_COMM_WORLD," %12f %12g\n",re,error);
134: }
135: }
136: PetscPrintf(PETSC_COMM_WORLD,"\n" );
137: }
138:
139: /*
140: Free work space
141: */
142: EPSDestroy(eps);
143: MatDestroy(A);
144: SlepcFinalize();
145: return 0;
146: }
148: /*
149: Compute the matrix vector multiplication y<---T*x where T is a nx by nx
150: tridiagonal matrix with DD on the diagonal, DL on the subdiagonal, and
151: DU on the superdiagonal.
152: */
153: static void tv( int nx, PetscScalar *x, PetscScalar *y )
154: {
155: PetscScalar dd, dl, du;
156: int j;
158: dd = 4.0;
159: dl = -1.0;
160: du = -1.0;
162: y[0] = dd*x[0] + du*x[1];
163: for( j=1; j<nx-1; j++ )
164: y[j] = dl*x[j-1] + dd*x[j] + du*x[j+1];
165: y[nx-1] = dl*x[nx-2] + dd*x[nx-1];
166: }
170: /*
171: Matrix-vector product subroutine for the 2D Laplacian.
173: The matrix used is the 2 dimensional discrete Laplacian on unit square with
174: zero Dirichlet boundary condition.
175:
176: Computes y <-- A*x, where A is the block tridiagonal matrix
177:
178: | T -I |
179: |-I T -I |
180: A = | -I T |
181: | ... -I|
182: | -I T|
183:
184: The subroutine TV is called to compute y<--T*x.
185: */
186: PetscErrorCode MatLaplacian2D_Mult( Mat A, Vec x, Vec y )
187: {
188: void *ctx;
190: int nx, lo, j, one=1;
191: PetscScalar *px, *py, dmone=-1.0;
192:
194: MatShellGetContext( A, &ctx );
195: nx = *(int *)ctx;
196: VecGetArray( x, &px );
197: VecGetArray( y, &py );
199: tv( nx, &px[0], &py[0] );
200: BLASaxpy_( &nx, &dmone, &px[nx], &one, &py[0], &one );
202: for( j=2; j<nx; j++ ) {
203: lo = (j-1)*nx;
204: tv( nx, &px[lo], &py[lo]);
205: BLASaxpy_( &nx, &dmone, &px[lo-nx], &one, &py[lo], &one );
206: BLASaxpy_( &nx, &dmone, &px[lo+nx], &one, &py[lo], &one );
207: }
209: lo = (nx-1)*nx;
210: tv( nx, &px[lo], &py[lo]);
211: BLASaxpy_( &nx, &dmone, &px[lo-nx], &one, &py[lo], &one );
213: VecRestoreArray( x, &px );
214: VecRestoreArray( y, &py );
215: return(0);
216: }