Actual source code: test8.c
slepc-3.17.0 2022-03-31
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[] = "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 MatMult_Laplacian2D(Mat A,Vec x,Vec y);
23: PetscErrorCode MatGetDiagonal_Laplacian2D(Mat A,Vec diag);
25: int main(int argc,char **argv)
26: {
27: Mat A; /* operator matrix */
28: EPS eps; /* eigenproblem solver context */
29: PetscReal tol=1000*PETSC_MACHINE_EPSILON;
30: PetscMPIInt size;
31: PetscInt N,n=10,nev;
33: SlepcInitialize(&argc,&argv,(char*)0,help);
34: MPI_Comm_size(PETSC_COMM_WORLD,&size);
37: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
38: N = n*n;
39: PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem (matrix-free version), N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,n);
41: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
42: Compute the operator matrix that defines the eigensystem, Ax=kx
43: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
45: MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,&n,&A);
46: MatShellSetOperation(A,MATOP_MULT,(void(*)(void))MatMult_Laplacian2D);
47: MatShellSetOperation(A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMult_Laplacian2D);
48: MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Laplacian2D);
50: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51: Create the eigensolver and set various options
52: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
54: /*
55: Create eigensolver context
56: */
57: EPSCreate(PETSC_COMM_WORLD,&eps);
59: /*
60: Set operators. In this case, it is a standard eigenvalue problem
61: */
62: EPSSetOperators(eps,A,NULL);
63: EPSSetProblemType(eps,EPS_HEP);
64: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
66: /*
67: Set solver parameters at runtime
68: */
69: EPSSetFromOptions(eps);
71: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72: Solve the eigensystem
73: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
75: EPSSolve(eps);
76: EPSGetDimensions(eps,&nev,NULL,NULL);
77: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
79: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80: Display solution and clean up
81: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
84: EPSDestroy(&eps);
85: MatDestroy(&A);
86: SlepcFinalize();
87: return 0;
88: }
90: /*
91: Compute the matrix vector multiplication y<---T*x where T is a nx by nx
92: tridiagonal matrix with DD on the diagonal, DL on the subdiagonal, and
93: DU on the superdiagonal.
94: */
95: static void tv(int nx,const PetscScalar *x,PetscScalar *y)
96: {
97: PetscScalar dd,dl,du;
98: int j;
100: dd = 4.0;
101: dl = -1.0;
102: du = -1.0;
104: y[0] = dd*x[0] + du*x[1];
105: for (j=1;j<nx-1;j++)
106: y[j] = dl*x[j-1] + dd*x[j] + du*x[j+1];
107: y[nx-1] = dl*x[nx-2] + dd*x[nx-1];
108: }
110: /*
111: Matrix-vector product subroutine for the 2D Laplacian.
113: The matrix used is the 2 dimensional discrete Laplacian on unit square with
114: zero Dirichlet boundary condition.
116: Computes y <-- A*x, where A is the block tridiagonal matrix
118: | T -I |
119: |-I T -I |
120: A = | -I T |
121: | ... -I|
122: | -I T|
124: The subroutine TV is called to compute y<--T*x.
125: */
126: PetscErrorCode MatMult_Laplacian2D(Mat A,Vec x,Vec y)
127: {
128: void *ctx;
129: int nx,lo,i,j;
130: const PetscScalar *px;
131: PetscScalar *py;
134: MatShellGetContext(A,&ctx);
135: nx = *(int*)ctx;
136: VecGetArrayRead(x,&px);
137: VecGetArray(y,&py);
139: tv(nx,&px[0],&py[0]);
140: for (i=0;i<nx;i++) py[i] -= px[nx+i];
142: for (j=2;j<nx;j++) {
143: lo = (j-1)*nx;
144: tv(nx,&px[lo],&py[lo]);
145: for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i] + px[lo+nx+i];
146: }
148: lo = (nx-1)*nx;
149: tv(nx,&px[lo],&py[lo]);
150: for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i];
152: VecRestoreArrayRead(x,&px);
153: VecRestoreArray(y,&py);
154: PetscFunctionReturn(0);
155: }
157: PetscErrorCode MatGetDiagonal_Laplacian2D(Mat A,Vec diag)
158: {
160: VecSet(diag,4.0);
161: PetscFunctionReturn(0);
162: }
164: /*TEST
166: testset:
167: args: -n 20 -eps_nev 4 -eps_ncv 11 -eps_max_it 40000
168: requires: !single
169: output_file: output/test8_1.out
170: test:
171: suffix: 1
172: args: -eps_type {{power subspace arnoldi}}
173: test:
174: suffix: 1_lanczos
175: args: -eps_type lanczos -eps_lanczos_reorthog local
176: test:
177: suffix: 1_lapack
178: args: -eps_type lapack
179: timeoutfactor: 2
180: test:
181: suffix: 1_elemental
182: args: -eps_type elemental
183: requires: elemental
184: test:
185: suffix: 1_krylovschur_vecs
186: args: -bv_type vecs -bv_orthog_refine always -eps_ncv 10
187: test:
188: suffix: 1_jd
189: args: -eps_type jd -eps_jd_blocksize 3
190: test:
191: suffix: 1_gd
192: args: -eps_type gd -eps_gd_blocksize 3 -eps_tol 1e-8
193: test:
194: suffix: 1_gd2
195: args: -eps_type gd -eps_gd_double_expansion
196: test:
197: suffix: 1_primme
198: args: -eps_type primme -eps_conv_abs
199: requires: primme
201: testset:
202: args: -eps_nev 4 -eps_smallest_real -eps_max_it 600
203: output_file: output/test8_2.out
204: test:
205: suffix: 2
206: args: -eps_type {{rqcg lobpcg}}
207: test:
208: suffix: 2_lanczos
209: args: -eps_type lanczos -eps_lanczos_reorthog local
210: test:
211: suffix: 2_arpack
212: args: -eps_type arpack -eps_ncv 6
213: requires: arpack !single
214: test:
215: suffix: 2_primme
216: args: -eps_type primme -eps_conv_abs -eps_primme_method lobpcg_orthobasisw -eps_ncv 24
217: requires: primme
219: testset:
220: args: -eps_nev 12 -eps_mpd 9 -eps_smallest_real -eps_max_it 1000
221: output_file: output/test8_3.out
222: test:
223: suffix: 3_rqcg
224: args: -eps_type rqcg
225: test:
226: suffix: 3_lanczos
227: args: -eps_type lanczos -eps_lanczos_reorthog local
228: test:
229: suffix: 3_lobpcg
230: args: -eps_type lobpcg -eps_lobpcg_blocksize 3 -eps_lobpcg_locking 0 -st_ksp_type preonly -st_pc_type jacobi
231: requires: !__float128
232: test:
233: suffix: 3_lobpcg_quad
234: args: -eps_type lobpcg -eps_lobpcg_blocksize 3 -eps_lobpcg_locking 0 -st_ksp_type preonly -st_pc_type jacobi -eps_tol 1e-25
235: requires: __float128
236: TEST*/