Actual source code: ex28.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[] = "A quadratic eigenproblem defined using shell matrices.\n\n"
23: "The command line options are:\n"
24: " -n <n>, where <n> = number of grid subdivisions in x and y dimensions.\n\n";
26: #include <slepcpep.h>
28: /*
29: User-defined routines
30: */
31: PetscErrorCode MatMult_Laplacian2D(Mat A,Vec x,Vec y);
32: PetscErrorCode MatGetDiagonal_Laplacian2D(Mat A,Vec diag);
33: PetscErrorCode MatMult_Zero(Mat A,Vec x,Vec y);
34: PetscErrorCode MatGetDiagonal_Zero(Mat A,Vec diag);
35: PetscErrorCode MatMult_Identity(Mat A,Vec x,Vec y);
36: PetscErrorCode MatGetDiagonal_Identity(Mat A,Vec diag);
40: int main(int argc,char **argv)
41: {
42: Mat M,C,K,A[3]; /* problem matrices */
43: PEP pep; /* polynomial eigenproblem solver context */
44: PEPType type;
45: PetscInt N,n=10,nev;
46: PetscMPIInt size;
47: PetscBool terse;
49: ST st;
51: SlepcInitialize(&argc,&argv,(char*)0,help);
52: MPI_Comm_size(PETSC_COMM_WORLD,&size);
53: if (size != 1) SETERRQ(PETSC_COMM_WORLD,1,"This is a uniprocessor example only");
55: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
56: N = n*n;
57: PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem with shell matrices, N=%D (%Dx%D grid)\n\n",N,n,n);
59: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
61: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63: /* K is the 2-D Laplacian */
64: MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,&n,&K);
65: MatSetFromOptions(K);
66: MatShellSetOperation(K,MATOP_MULT,(void(*)())MatMult_Laplacian2D);
67: MatShellSetOperation(K,MATOP_MULT_TRANSPOSE,(void(*)())MatMult_Laplacian2D);
68: MatShellSetOperation(K,MATOP_GET_DIAGONAL,(void(*)())MatGetDiagonal_Laplacian2D);
70: /* C is the zero matrix */
71: MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,NULL,&C);
72: MatSetFromOptions(C);
73: MatShellSetOperation(C,MATOP_MULT,(void(*)())MatMult_Zero);
74: MatShellSetOperation(C,MATOP_MULT_TRANSPOSE,(void(*)())MatMult_Zero);
75: MatShellSetOperation(C,MATOP_GET_DIAGONAL,(void(*)())MatGetDiagonal_Zero);
77: /* M is the identity matrix */
78: MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,NULL,&M);
79: MatSetFromOptions(M);
80: MatShellSetOperation(M,MATOP_MULT,(void(*)())MatMult_Identity);
81: MatShellSetOperation(M,MATOP_MULT_TRANSPOSE,(void(*)())MatMult_Identity);
82: MatShellSetOperation(M,MATOP_GET_DIAGONAL,(void(*)())MatGetDiagonal_Identity);
84: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85: Create the eigensolver and set various options
86: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
88: /*
89: Create eigensolver context
90: */
91: PEPCreate(PETSC_COMM_WORLD,&pep);
93: /*
94: Set matrices and problem type
95: */
96: A[0] = K; A[1] = C; A[2] = M;
97: PEPSetOperators(pep,3,A);
98: PEPGetST(pep,&st);
99: STSetMatMode(st,ST_MATMODE_SHELL);
101: /*
102: Set solver parameters at runtime
103: */
104: PEPSetFromOptions(pep);
106: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107: Solve the eigensystem
108: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110: PEPSolve(pep);
112: /*
113: Optional: Get some information from the solver and display it
114: */
115: PEPGetType(pep,&type);
116: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
117: PEPGetDimensions(pep,&nev,NULL,NULL);
118: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
120: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121: Display solution and clean up
122: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124: /* show detailed info unless -terse option is given by user */
125: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
126: if (terse) {
127: PEPErrorView(pep,PEP_ERROR_RELATIVE,NULL);
128: } else {
129: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
130: PEPReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
131: PEPErrorView(pep,PEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
132: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
133: }
134: PEPDestroy(&pep);
135: MatDestroy(&M);
136: MatDestroy(&C);
137: MatDestroy(&K);
138: SlepcFinalize();
139: return ierr;
140: }
142: /*
143: Compute the matrix vector multiplication y<---T*x where T is a nx by nx
144: tridiagonal matrix with DD on the diagonal, DL on the subdiagonal, and
145: DU on the superdiagonal.
146: */
147: static void tv(int nx,const PetscScalar *x,PetscScalar *y)
148: {
149: PetscScalar dd,dl,du;
150: int j;
152: dd = 4.0;
153: dl = -1.0;
154: du = -1.0;
156: y[0] = dd*x[0] + du*x[1];
157: for (j=1;j<nx-1;j++)
158: y[j] = dl*x[j-1] + dd*x[j] + du*x[j+1];
159: y[nx-1] = dl*x[nx-2] + dd*x[nx-1];
160: }
164: /*
165: Matrix-vector product subroutine for the 2D Laplacian.
167: The matrix used is the 2 dimensional discrete Laplacian on unit square with
168: zero Dirichlet boundary condition.
170: Computes y <-- A*x, where A is the block tridiagonal matrix
172: | T -I |
173: |-I T -I |
174: A = | -I T |
175: | ... -I|
176: | -I T|
178: The subroutine TV is called to compute y<--T*x.
179: */
180: PetscErrorCode MatMult_Laplacian2D(Mat A,Vec x,Vec y)
181: {
182: void *ctx;
183: int nx,lo,i,j;
184: const PetscScalar *px;
185: PetscScalar *py;
186: PetscErrorCode ierr;
189: MatShellGetContext(A,&ctx);
190: nx = *(int*)ctx;
191: VecGetArrayRead(x,&px);
192: VecGetArray(y,&py);
194: tv(nx,&px[0],&py[0]);
195: for (i=0;i<nx;i++) py[i] -= px[nx+i];
197: for (j=2;j<nx;j++) {
198: lo = (j-1)*nx;
199: tv(nx,&px[lo],&py[lo]);
200: for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i] + px[lo+nx+i];
201: }
203: lo = (nx-1)*nx;
204: tv(nx,&px[lo],&py[lo]);
205: for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i];
207: VecRestoreArrayRead(x,&px);
208: VecRestoreArray(y,&py);
209: return(0);
210: }
214: PetscErrorCode MatGetDiagonal_Laplacian2D(Mat A,Vec diag)
215: {
219: VecSet(diag,4.0);
220: return(0);
221: }
225: /*
226: Matrix-vector product subroutine for the Null matrix.
227: */
228: PetscErrorCode MatMult_Zero(Mat A,Vec x,Vec y)
229: {
233: VecSet(y,0.0);
234: return(0);
235: }
239: PetscErrorCode MatGetDiagonal_Zero(Mat A,Vec diag)
240: {
244: VecSet(diag,0.0);
245: return(0);
246: }
250: /*
251: Matrix-vector product subroutine for the Identity matrix.
252: */
253: PetscErrorCode MatMult_Identity(Mat A,Vec x,Vec y)
254: {
255: PetscErrorCode ierr;
258: VecCopy(x,y);
259: return(0);
260: }
264: PetscErrorCode MatGetDiagonal_Identity(Mat A,Vec diag)
265: {
269: VecSet(diag,1.0);
270: return(0);
271: }