Actual source code: ex28.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: static char help[] = "A quadratic eigenproblem defined using shell matrices.\n\n"
12: "The command line options are:\n"
13: " -n <n>, where <n> = number of grid subdivisions in x and y dimensions.\n\n";
15: #include <slepcpep.h>
17: /*
18: User-defined routines
19: */
20: PetscErrorCode MatMult_Laplacian2D(Mat A,Vec x,Vec y);
21: PetscErrorCode MatGetDiagonal_Laplacian2D(Mat A,Vec diag);
22: PetscErrorCode MatMult_Zero(Mat A,Vec x,Vec y);
23: PetscErrorCode MatGetDiagonal_Zero(Mat A,Vec diag);
24: PetscErrorCode MatMult_Identity(Mat A,Vec x,Vec y);
25: PetscErrorCode MatGetDiagonal_Identity(Mat A,Vec diag);
27: int main(int argc,char **argv)
28: {
29: Mat M,C,K,A[3]; /* problem matrices */
30: PEP pep; /* polynomial eigenproblem solver context */
31: PEPType type;
32: PetscInt N,n=10,nev;
33: PetscMPIInt size;
34: PetscBool terse;
36: ST st;
38: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
39: MPI_Comm_size(PETSC_COMM_WORLD,&size);
40: if (size != 1) SETERRQ(PETSC_COMM_WORLD,1,"This is a uniprocessor example only");
42: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
43: N = n*n;
44: PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem with shell matrices, N=%D (%Dx%D grid)\n\n",N,n,n);
46: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
48: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50: /* K is the 2-D Laplacian */
51: MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,&n,&K);
52: MatShellSetOperation(K,MATOP_MULT,(void(*)(void))MatMult_Laplacian2D);
53: MatShellSetOperation(K,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMult_Laplacian2D);
54: MatShellSetOperation(K,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Laplacian2D);
56: /* C is the zero matrix */
57: MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,NULL,&C);
58: MatShellSetOperation(C,MATOP_MULT,(void(*)(void))MatMult_Zero);
59: MatShellSetOperation(C,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMult_Zero);
60: MatShellSetOperation(C,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Zero);
62: /* M is the identity matrix */
63: MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,NULL,&M);
64: MatShellSetOperation(M,MATOP_MULT,(void(*)(void))MatMult_Identity);
65: MatShellSetOperation(M,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMult_Identity);
66: MatShellSetOperation(M,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Identity);
68: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69: Create the eigensolver and set various options
70: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
72: /*
73: Create eigensolver context
74: */
75: PEPCreate(PETSC_COMM_WORLD,&pep);
77: /*
78: Set matrices and problem type
79: */
80: A[0] = K; A[1] = C; A[2] = M;
81: PEPSetOperators(pep,3,A);
82: PEPGetST(pep,&st);
83: STSetMatMode(st,ST_MATMODE_SHELL);
85: /*
86: Set solver parameters at runtime
87: */
88: PEPSetFromOptions(pep);
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Solve the eigensystem
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94: PEPSolve(pep);
96: /*
97: Optional: Get some information from the solver and display it
98: */
99: PEPGetType(pep,&type);
100: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
101: PEPGetDimensions(pep,&nev,NULL,NULL);
102: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
104: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105: Display solution and clean up
106: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108: /* show detailed info unless -terse option is given by user */
109: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
110: if (terse) {
111: PEPErrorView(pep,PEP_ERROR_RELATIVE,NULL);
112: } else {
113: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
114: PEPReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
115: PEPErrorView(pep,PEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
116: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
117: }
118: PEPDestroy(&pep);
119: MatDestroy(&M);
120: MatDestroy(&C);
121: MatDestroy(&K);
122: SlepcFinalize();
123: return ierr;
124: }
126: /*
127: Compute the matrix vector multiplication y<---T*x where T is a nx by nx
128: tridiagonal matrix with DD on the diagonal, DL on the subdiagonal, and
129: DU on the superdiagonal.
130: */
131: static void tv(int nx,const PetscScalar *x,PetscScalar *y)
132: {
133: PetscScalar dd,dl,du;
134: int j;
136: dd = 4.0;
137: dl = -1.0;
138: du = -1.0;
140: y[0] = dd*x[0] + du*x[1];
141: for (j=1;j<nx-1;j++)
142: y[j] = dl*x[j-1] + dd*x[j] + du*x[j+1];
143: y[nx-1] = dl*x[nx-2] + dd*x[nx-1];
144: }
146: /*
147: Matrix-vector product subroutine for the 2D Laplacian.
149: The matrix used is the 2 dimensional discrete Laplacian on unit square with
150: zero Dirichlet boundary condition.
152: Computes y <-- A*x, where A is the block tridiagonal matrix
154: | T -I |
155: |-I T -I |
156: A = | -I T |
157: | ... -I|
158: | -I T|
160: The subroutine TV is called to compute y<--T*x.
161: */
162: PetscErrorCode MatMult_Laplacian2D(Mat A,Vec x,Vec y)
163: {
164: void *ctx;
165: int nx,lo,i,j;
166: const PetscScalar *px;
167: PetscScalar *py;
168: PetscErrorCode ierr;
171: MatShellGetContext(A,&ctx);
172: nx = *(int*)ctx;
173: VecGetArrayRead(x,&px);
174: VecGetArray(y,&py);
176: tv(nx,&px[0],&py[0]);
177: for (i=0;i<nx;i++) py[i] -= px[nx+i];
179: for (j=2;j<nx;j++) {
180: lo = (j-1)*nx;
181: tv(nx,&px[lo],&py[lo]);
182: for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i] + px[lo+nx+i];
183: }
185: lo = (nx-1)*nx;
186: tv(nx,&px[lo],&py[lo]);
187: for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i];
189: VecRestoreArrayRead(x,&px);
190: VecRestoreArray(y,&py);
191: return(0);
192: }
194: PetscErrorCode MatGetDiagonal_Laplacian2D(Mat A,Vec diag)
195: {
199: VecSet(diag,4.0);
200: return(0);
201: }
203: /*
204: Matrix-vector product subroutine for the Null matrix.
205: */
206: PetscErrorCode MatMult_Zero(Mat A,Vec x,Vec y)
207: {
211: VecSet(y,0.0);
212: return(0);
213: }
215: PetscErrorCode MatGetDiagonal_Zero(Mat A,Vec diag)
216: {
220: VecSet(diag,0.0);
221: return(0);
222: }
224: /*
225: Matrix-vector product subroutine for the Identity matrix.
226: */
227: PetscErrorCode MatMult_Identity(Mat A,Vec x,Vec y)
228: {
229: PetscErrorCode ierr;
232: VecCopy(x,y);
233: return(0);
234: }
236: PetscErrorCode MatGetDiagonal_Identity(Mat A,Vec diag)
237: {
241: VecSet(diag,1.0);
242: return(0);
243: }
245: /*TEST
247: test:
248: suffix: 1
249: args: -pep_type {{toar qarnoldi linear}} -pep_nev 4 -terse
250: filter: grep -v Solution
251: filter_output: grep -v Solution
252: requires: !complex !single
254: TEST*/