Actual source code: ex28.c
slepc-3.18.3 2023-03-24
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;
38: SlepcInitialize(&argc,&argv,(char*)0,help);
39: MPI_Comm_size(PETSC_COMM_WORLD,&size);
42: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
43: N = n*n;
44: PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem with shell matrices, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " 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: %" PetscInt_FMT "\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) PEPErrorView(pep,PEP_ERROR_RELATIVE,NULL);
111: else {
112: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
113: PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
114: PEPErrorView(pep,PEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
115: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
116: }
117: PEPDestroy(&pep);
118: MatDestroy(&M);
119: MatDestroy(&C);
120: MatDestroy(&K);
121: SlepcFinalize();
122: return 0;
123: }
125: /*
126: Compute the matrix vector multiplication y<---T*x where T is a nx by nx
127: tridiagonal matrix with DD on the diagonal, DL on the subdiagonal, and
128: DU on the superdiagonal.
129: */
130: static void tv(int nx,const PetscScalar *x,PetscScalar *y)
131: {
132: PetscScalar dd,dl,du;
133: int j;
135: dd = 4.0;
136: dl = -1.0;
137: du = -1.0;
139: y[0] = dd*x[0] + du*x[1];
140: for (j=1;j<nx-1;j++)
141: y[j] = dl*x[j-1] + dd*x[j] + du*x[j+1];
142: y[nx-1] = dl*x[nx-2] + dd*x[nx-1];
143: }
145: /*
146: Matrix-vector product subroutine for the 2D Laplacian.
148: The matrix used is the 2 dimensional discrete Laplacian on unit square with
149: zero Dirichlet boundary condition.
151: Computes y <-- A*x, where A is the block tridiagonal matrix
153: | T -I |
154: |-I T -I |
155: A = | -I T |
156: | ... -I|
157: | -I T|
159: The subroutine TV is called to compute y<--T*x.
160: */
161: PetscErrorCode MatMult_Laplacian2D(Mat A,Vec x,Vec y)
162: {
163: void *ctx;
164: int nx,lo,i,j;
165: const PetscScalar *px;
166: PetscScalar *py;
169: MatShellGetContext(A,&ctx);
170: nx = *(int*)ctx;
171: VecGetArrayRead(x,&px);
172: VecGetArray(y,&py);
174: tv(nx,&px[0],&py[0]);
175: for (i=0;i<nx;i++) py[i] -= px[nx+i];
177: for (j=2;j<nx;j++) {
178: lo = (j-1)*nx;
179: tv(nx,&px[lo],&py[lo]);
180: for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i] + px[lo+nx+i];
181: }
183: lo = (nx-1)*nx;
184: tv(nx,&px[lo],&py[lo]);
185: for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i];
187: VecRestoreArrayRead(x,&px);
188: VecRestoreArray(y,&py);
189: return 0;
190: }
192: PetscErrorCode MatGetDiagonal_Laplacian2D(Mat A,Vec diag)
193: {
195: VecSet(diag,4.0);
196: return 0;
197: }
199: /*
200: Matrix-vector product subroutine for the Null matrix.
201: */
202: PetscErrorCode MatMult_Zero(Mat A,Vec x,Vec y)
203: {
205: VecSet(y,0.0);
206: return 0;
207: }
209: PetscErrorCode MatGetDiagonal_Zero(Mat A,Vec diag)
210: {
212: VecSet(diag,0.0);
213: return 0;
214: }
216: /*
217: Matrix-vector product subroutine for the Identity matrix.
218: */
219: PetscErrorCode MatMult_Identity(Mat A,Vec x,Vec y)
220: {
222: VecCopy(x,y);
223: return 0;
224: }
226: PetscErrorCode MatGetDiagonal_Identity(Mat A,Vec diag)
227: {
229: VecSet(diag,1.0);
230: return 0;
231: }
233: /*TEST
235: test:
236: suffix: 1
237: args: -pep_type {{toar qarnoldi linear}} -pep_nev 4 -terse
238: 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/"
240: TEST*/