Actual source code: test8.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[] = "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;
34: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
35: MPI_Comm_size(PETSC_COMM_WORLD,&size);
36: if (size != 1) SETERRQ(PETSC_COMM_WORLD,1,"This is a uniprocessor example only");
38: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
39: N = n*n;
40: PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem (matrix-free version), N=%D (%Dx%D grid)\n\n",N,n,n);
42: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43: Compute the operator matrix that defines the eigensystem, Ax=kx
44: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
46: MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,&n,&A);
47: MatShellSetOperation(A,MATOP_MULT,(void(*)(void))MatMult_Laplacian2D);
48: MatShellSetOperation(A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMult_Laplacian2D);
49: MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Laplacian2D);
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Create the eigensolver and set various options
53: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: /*
56: Create eigensolver context
57: */
58: EPSCreate(PETSC_COMM_WORLD,&eps);
60: /*
61: Set operators. In this case, it is a standard eigenvalue problem
62: */
63: EPSSetOperators(eps,A,NULL);
64: EPSSetProblemType(eps,EPS_HEP);
65: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
67: /*
68: Set solver parameters at runtime
69: */
70: EPSSetFromOptions(eps);
72: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73: Solve the eigensystem
74: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
76: EPSSolve(eps);
77: EPSGetDimensions(eps,&nev,NULL,NULL);
78: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
80: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81: Display solution and clean up
82: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
84: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
85: EPSDestroy(&eps);
86: MatDestroy(&A);
87: SlepcFinalize();
88: return ierr;
89: }
91: /*
92: Compute the matrix vector multiplication y<---T*x where T is a nx by nx
93: tridiagonal matrix with DD on the diagonal, DL on the subdiagonal, and
94: DU on the superdiagonal.
95: */
96: static void tv(int nx,const PetscScalar *x,PetscScalar *y)
97: {
98: PetscScalar dd,dl,du;
99: int j;
101: dd = 4.0;
102: dl = -1.0;
103: du = -1.0;
105: y[0] = dd*x[0] + du*x[1];
106: for (j=1;j<nx-1;j++)
107: y[j] = dl*x[j-1] + dd*x[j] + du*x[j+1];
108: y[nx-1] = dl*x[nx-2] + dd*x[nx-1];
109: }
111: /*
112: Matrix-vector product subroutine for the 2D Laplacian.
114: The matrix used is the 2 dimensional discrete Laplacian on unit square with
115: zero Dirichlet boundary condition.
117: Computes y <-- A*x, where A is the block tridiagonal matrix
119: | T -I |
120: |-I T -I |
121: A = | -I T |
122: | ... -I|
123: | -I T|
125: The subroutine TV is called to compute y<--T*x.
126: */
127: PetscErrorCode MatMult_Laplacian2D(Mat A,Vec x,Vec y)
128: {
129: void *ctx;
130: int nx,lo,i,j;
131: const PetscScalar *px;
132: PetscScalar *py;
133: PetscErrorCode ierr;
136: MatShellGetContext(A,&ctx);
137: nx = *(int*)ctx;
138: VecGetArrayRead(x,&px);
139: VecGetArray(y,&py);
141: tv(nx,&px[0],&py[0]);
142: for (i=0;i<nx;i++) py[i] -= px[nx+i];
144: for (j=2;j<nx;j++) {
145: lo = (j-1)*nx;
146: tv(nx,&px[lo],&py[lo]);
147: for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i] + px[lo+nx+i];
148: }
150: lo = (nx-1)*nx;
151: tv(nx,&px[lo],&py[lo]);
152: for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i];
154: VecRestoreArrayRead(x,&px);
155: VecRestoreArray(y,&py);
156: return(0);
157: }
159: PetscErrorCode MatGetDiagonal_Laplacian2D(Mat A,Vec diag)
160: {
164: VecSet(diag,4.0);
165: return(0);
166: }
168: /*TEST
170: testset:
171: args: -n 20 -eps_nev 4 -eps_ncv 11 -eps_max_it 40000
172: requires: !single
173: output_file: output/test8_1.out
174: test:
175: suffix: 1
176: args: -eps_type {{krylovschur power subspace arnoldi lanczos lapack}}
177: test:
178: suffix: 1_krylovschur_vecs
179: args: -bv_type vecs -bv_orthog_refine always -eps_ncv 12
180: test:
181: suffix: 1_jd
182: args: -eps_type jd -eps_jd_blocksize 3
183: test:
184: suffix: 1_gd
185: args: -eps_type gd -eps_gd_blocksize 3 -eps_tol 1e-8
186: test:
187: suffix: 1_gd2
188: args: -eps_type gd -eps_gd_double_expansion
189: test:
190: suffix: 1_primme
191: args: -eps_type primme -eps_conv_abs
192: requires: primme
194: testset:
195: args: -eps_nev 4 -eps_smallest_real -eps_max_it 500
196: output_file: output/test8_2.out
197: test:
198: suffix: 2
199: args: -eps_type {{rqcg lobpcg lanczos}}
200: requires: !single
201: test:
202: suffix: 2_single
203: args: -eps_type {{rqcg lobpcg lanczos}} -eps_tol 1e-5
204: requires: single
205: test:
206: suffix: 2_arpack
207: args: -eps_type arpack -eps_ncv 6
208: requires: arpack !single
209: test:
210: suffix: 2_blzpack
211: args: -eps_type blzpack
212: requires: blzpack
213: test:
214: suffix: 2_blopex
215: args: -eps_type blopex
216: requires: blopex
218: testset:
219: args: -eps_nev 12 -eps_mpd 9 -eps_smallest_real -eps_max_it 1000
220: output_file: output/test8_3.out
221: test:
222: suffix: 3
223: args: -eps_type {{rqcg lanczos}}
224: requires: double
225: test:
226: suffix: 3_lobpcg
227: args: -eps_type lobpcg -eps_lobpcg_blocksize 3 -eps_lobpcg_locking 0 -st_ksp_type preonly -st_pc_type jacobi
228: requires: double
229: test:
230: suffix: 3_single
231: args: -eps_type {{rqcg lanczos}} -eps_tol 1e-5
232: requires: single
233: test:
234: suffix: 3_lobpcg_single
235: args: -eps_type lobpcg -eps_lobpcg_blocksize 3 -eps_lobpcg_locking 0 -st_ksp_type preonly -st_pc_type jacobi -eps_tol 1e-5
236: requires: single
237: test:
238: suffix: 3_quad
239: args: -eps_type {{rqcg lanczos}} -eps_tol 1e-25
240: requires: __float128
241: test:
242: suffix: 3_lobpcg_quad
243: args: -eps_type lobpcg -eps_lobpcg_blocksize 3 -eps_lobpcg_locking 0 -st_ksp_type preonly -st_pc_type jacobi -eps_tol 1e-25
244: requires: __float128
245: TEST*/