Actual source code: ex28.c
slepc-3.17.1 2022-04-11
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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;
35: ST st;
37: SlepcInitialize(&argc,&argv,(char*)0,help);
38: MPI_Comm_size(PETSC_COMM_WORLD,&size);
41: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
42: N = n*n;
43: PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem with shell matrices, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,n);
45: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
46: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
47: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
49: /* K is the 2-D Laplacian */
50: MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,&n,&K);
51: MatShellSetOperation(K,MATOP_MULT,(void(*)(void))MatMult_Laplacian2D);
52: MatShellSetOperation(K,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMult_Laplacian2D);
53: MatShellSetOperation(K,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Laplacian2D);
55: /* C is the zero matrix */
56: MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,NULL,&C);
57: MatShellSetOperation(C,MATOP_MULT,(void(*)(void))MatMult_Zero);
58: MatShellSetOperation(C,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMult_Zero);
59: MatShellSetOperation(C,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Zero);
61: /* M is the identity matrix */
62: MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,NULL,&M);
63: MatShellSetOperation(M,MATOP_MULT,(void(*)(void))MatMult_Identity);
64: MatShellSetOperation(M,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMult_Identity);
65: MatShellSetOperation(M,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Identity);
67: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
68: Create the eigensolver and set various options
69: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71: /*
72: Create eigensolver context
73: */
74: PEPCreate(PETSC_COMM_WORLD,&pep);
76: /*
77: Set matrices and problem type
78: */
79: A[0] = K; A[1] = C; A[2] = M;
80: PEPSetOperators(pep,3,A);
81: PEPGetST(pep,&st);
82: STSetMatMode(st,ST_MATMODE_SHELL);
84: /*
85: Set solver parameters at runtime
86: */
87: PEPSetFromOptions(pep);
89: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
90: Solve the eigensystem
91: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93: PEPSolve(pep);
95: /*
96: Optional: Get some information from the solver and display it
97: */
98: PEPGetType(pep,&type);
99: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
100: PEPGetDimensions(pep,&nev,NULL,NULL);
101: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
103: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104: Display solution and clean up
105: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107: /* show detailed info unless -terse option is given by user */
108: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
109: if (terse) PEPErrorView(pep,PEP_ERROR_RELATIVE,NULL);
110: else {
111: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
112: PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
113: PEPErrorView(pep,PEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
114: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
115: }
116: PEPDestroy(&pep);
117: MatDestroy(&M);
118: MatDestroy(&C);
119: MatDestroy(&K);
120: SlepcFinalize();
121: return 0;
122: }
124: /*
125: Compute the matrix vector multiplication y<---T*x where T is a nx by nx
126: tridiagonal matrix with DD on the diagonal, DL on the subdiagonal, and
127: DU on the superdiagonal.
128: */
129: static void tv(int nx,const PetscScalar *x,PetscScalar *y)
130: {
131: PetscScalar dd,dl,du;
132: int j;
134: dd = 4.0;
135: dl = -1.0;
136: du = -1.0;
138: y[0] = dd*x[0] + du*x[1];
139: for (j=1;j<nx-1;j++)
140: y[j] = dl*x[j-1] + dd*x[j] + du*x[j+1];
141: y[nx-1] = dl*x[nx-2] + dd*x[nx-1];
142: }
144: /*
145: Matrix-vector product subroutine for the 2D Laplacian.
147: The matrix used is the 2 dimensional discrete Laplacian on unit square with
148: zero Dirichlet boundary condition.
150: Computes y <-- A*x, where A is the block tridiagonal matrix
152: | T -I |
153: |-I T -I |
154: A = | -I T |
155: | ... -I|
156: | -I T|
158: The subroutine TV is called to compute y<--T*x.
159: */
160: PetscErrorCode MatMult_Laplacian2D(Mat A,Vec x,Vec y)
161: {
162: void *ctx;
163: int nx,lo,i,j;
164: const PetscScalar *px;
165: PetscScalar *py;
168: MatShellGetContext(A,&ctx);
169: nx = *(int*)ctx;
170: VecGetArrayRead(x,&px);
171: VecGetArray(y,&py);
173: tv(nx,&px[0],&py[0]);
174: for (i=0;i<nx;i++) py[i] -= px[nx+i];
176: for (j=2;j<nx;j++) {
177: lo = (j-1)*nx;
178: tv(nx,&px[lo],&py[lo]);
179: for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i] + px[lo+nx+i];
180: }
182: lo = (nx-1)*nx;
183: tv(nx,&px[lo],&py[lo]);
184: for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i];
186: VecRestoreArrayRead(x,&px);
187: VecRestoreArray(y,&py);
188: PetscFunctionReturn(0);
189: }
191: PetscErrorCode MatGetDiagonal_Laplacian2D(Mat A,Vec diag)
192: {
194: VecSet(diag,4.0);
195: PetscFunctionReturn(0);
196: }
198: /*
199: Matrix-vector product subroutine for the Null matrix.
200: */
201: PetscErrorCode MatMult_Zero(Mat A,Vec x,Vec y)
202: {
204: VecSet(y,0.0);
205: PetscFunctionReturn(0);
206: }
208: PetscErrorCode MatGetDiagonal_Zero(Mat A,Vec diag)
209: {
211: VecSet(diag,0.0);
212: PetscFunctionReturn(0);
213: }
215: /*
216: Matrix-vector product subroutine for the Identity matrix.
217: */
218: PetscErrorCode MatMult_Identity(Mat A,Vec x,Vec y)
219: {
221: VecCopy(x,y);
222: PetscFunctionReturn(0);
223: }
225: PetscErrorCode MatGetDiagonal_Identity(Mat A,Vec diag)
226: {
228: VecSet(diag,1.0);
229: PetscFunctionReturn(0);
230: }
232: /*TEST
234: test:
235: suffix: 1
236: args: -pep_type {{toar qarnoldi linear}} -pep_nev 4 -terse
237: filter: grep -v Solution | sed -e "s/2.7996[1-8]i/2.79964i/g" | sed -e "s/2.7570[5-9]i/2.75708i/g" | sed -e "s/0.00000-2.79964i, 0.00000+2.79964i/0.00000+2.79964i, 0.00000-2.79964i/" | sed -e "s/0.00000-2.75708i, 0.00000+2.75708i/0.00000+2.75708i, 0.00000-2.75708i/"
239: TEST*/