Line | Branch | Exec | Source |
---|---|---|---|
1 | /* | ||
2 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
3 | SLEPc - Scalable Library for Eigenvalue Problem Computations | ||
4 | Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain | ||
5 | |||
6 | This file is part of SLEPc. | ||
7 | SLEPc is distributed under a 2-clause BSD license (see LICENSE). | ||
8 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
9 | */ | ||
10 | |||
11 | static char help[] = "Solves the same eigenproblem as in example ex2, but using a shell matrix.\n\n" | ||
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"; | ||
15 | |||
16 | #include <slepceps.h> | ||
17 | #include <petscblaslapack.h> | ||
18 | |||
19 | /* | ||
20 | User-defined routines | ||
21 | */ | ||
22 | PetscErrorCode MatMult_Laplacian2D(Mat A,Vec x,Vec y); | ||
23 | PetscErrorCode MatGetDiagonal_Laplacian2D(Mat A,Vec diag); | ||
24 | |||
25 | 169 | int main(int argc,char **argv) | |
26 | { | ||
27 | 169 | Mat A; /* operator matrix */ | |
28 | 169 | EPS eps; /* eigenproblem solver context */ | |
29 | 169 | PetscReal tol=1000*PETSC_MACHINE_EPSILON; | |
30 | 169 | PetscMPIInt size; | |
31 | 169 | PetscInt N,n=10,nev; | |
32 | |||
33 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
169 | PetscFunctionBeginUser; |
34 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
35 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
169 | PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); |
36 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
169 | PetscCheck(size==1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only"); |
37 | |||
38 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
39 | 169 | N = n*n; | |
40 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(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 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
43 | Compute the operator matrix that defines the eigensystem, Ax=kx | ||
44 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
45 | |||
46 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,&n,&A)); |
47 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(MatShellSetOperation(A,MATOP_MULT,(PetscErrorCodeFn*)MatMult_Laplacian2D)); |
48 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(MatShellSetOperation(A,MATOP_MULT_TRANSPOSE,(PetscErrorCodeFn*)MatMult_Laplacian2D)); |
49 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(MatShellSetOperation(A,MATOP_GET_DIAGONAL,(PetscErrorCodeFn*)MatGetDiagonal_Laplacian2D)); |
50 | |||
51 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
52 | Create the eigensolver and set various options | ||
53 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
54 | |||
55 | /* | ||
56 | Create eigensolver context | ||
57 | */ | ||
58 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps)); |
59 | |||
60 | /* | ||
61 | Set operators. In this case, it is a standard eigenvalue problem | ||
62 | */ | ||
63 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(EPSSetOperators(eps,A,NULL)); |
64 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(EPSSetProblemType(eps,EPS_HEP)); |
65 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(EPSSetTolerances(eps,tol,PETSC_CURRENT)); |
66 | |||
67 | /* | ||
68 | Set solver parameters at runtime | ||
69 | */ | ||
70 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(EPSSetFromOptions(eps)); |
71 | |||
72 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
73 | Solve the eigensystem | ||
74 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
75 | |||
76 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(EPSSolve(eps)); |
77 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL)); |
78 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev)); |
79 | |||
80 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
81 | Display solution and clean up | ||
82 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
83 | |||
84 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL)); |
85 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(EPSDestroy(&eps)); |
86 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
169 | PetscCall(MatDestroy(&A)); |
87 |
3/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
169 | PetscCall(SlepcFinalize()); |
88 | return 0; | ||
89 | } | ||
90 | |||
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 | 9106310 | static void tv(int nx,const PetscScalar *x,PetscScalar *y) | |
97 | { | ||
98 | 9106310 | PetscScalar dd,dl,du; | |
99 | 9106310 | int j; | |
100 | |||
101 | 9106310 | dd = 4.0; | |
102 | 9106310 | dl = -1.0; | |
103 | 9106310 | du = -1.0; | |
104 | |||
105 | 9106310 | y[0] = dd*x[0] + du*x[1]; | |
106 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
150498190 | for (j=1;j<nx-1;j++) |
107 | 141391880 | y[j] = dl*x[j-1] + dd*x[j] + du*x[j+1]; | |
108 | 9106310 | y[nx-1] = dl*x[nx-2] + dd*x[nx-1]; | |
109 | 9106310 | } | |
110 | |||
111 | /* | ||
112 | Matrix-vector product subroutine for the 2D Laplacian. | ||
113 | |||
114 | The matrix used is the 2 dimensional discrete Laplacian on unit square with | ||
115 | zero Dirichlet boundary condition. | ||
116 | |||
117 | Computes y <-- A*x, where A is the block tridiagonal matrix | ||
118 | |||
119 | | T -I | | ||
120 | |-I T -I | | ||
121 | A = | -I T | | ||
122 | | ... -I| | ||
123 | | -I T| | ||
124 | |||
125 | The subroutine TV is called to compute y<--T*x. | ||
126 | */ | ||
127 | 567924 | PetscErrorCode MatMult_Laplacian2D(Mat A,Vec x,Vec y) | |
128 | { | ||
129 | 567924 | void *ctx; | |
130 | 567924 | int nx,lo,i,j; | |
131 | 567924 | const PetscScalar *px; | |
132 | 567924 | PetscScalar *py; | |
133 | |||
134 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
567924 | PetscFunctionBeginUser; |
135 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
567924 | PetscCall(MatShellGetContext(A,&ctx)); |
136 | 567924 | nx = *(int*)ctx; | |
137 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
567924 | PetscCall(VecGetArrayRead(x,&px)); |
138 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
567924 | PetscCall(VecGetArray(y,&py)); |
139 | |||
140 | 567924 | tv(nx,&px[0],&py[0]); | |
141 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
10242158 | for (i=0;i<nx;i++) py[i] -= px[nx+i]; |
142 | |||
143 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
8538386 | for (j=2;j<nx;j++) { |
144 | 7970462 | lo = (j-1)*nx; | |
145 | 7970462 | tv(nx,&px[lo],&py[lo]); | |
146 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
157332804 | for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i] + px[lo+nx+i]; |
147 | } | ||
148 | |||
149 | 567924 | lo = (nx-1)*nx; | |
150 | 567924 | tv(nx,&px[lo],&py[lo]); | |
151 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
10242158 | for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i]; |
152 | |||
153 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
567924 | PetscCall(VecRestoreArrayRead(x,&px)); |
154 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
567924 | PetscCall(VecRestoreArray(y,&py)); |
155 |
5/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
121648 | PetscFunctionReturn(PETSC_SUCCESS); |
156 | } | ||
157 | |||
158 | 10 | PetscErrorCode MatGetDiagonal_Laplacian2D(Mat A,Vec diag) | |
159 | { | ||
160 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBeginUser; |
161 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(VecSet(diag,4.0)); |
162 |
5/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
2 | PetscFunctionReturn(PETSC_SUCCESS); |
163 | } | ||
164 | |||
165 | /*TEST | ||
166 | |||
167 | testset: | ||
168 | args: -n 20 -eps_nev 4 -eps_ncv 11 -eps_max_it 40000 | ||
169 | requires: !single | ||
170 | output_file: output/test8_1.out | ||
171 | test: | ||
172 | suffix: 1 | ||
173 | args: -eps_type {{power subspace arnoldi}} | ||
174 | test: | ||
175 | suffix: 1_lanczos | ||
176 | args: -eps_type lanczos -eps_lanczos_reorthog local | ||
177 | test: | ||
178 | suffix: 1_lapack | ||
179 | args: -eps_type lapack | ||
180 | timeoutfactor: 2 | ||
181 | test: | ||
182 | suffix: 1_elemental | ||
183 | args: -eps_type elemental | ||
184 | requires: elemental | ||
185 | test: | ||
186 | suffix: 1_krylovschur_vecs | ||
187 | args: -bv_type vecs -bv_orthog_refine always -eps_ncv 10 -vec_mdot_use_gemv 0 | ||
188 | test: | ||
189 | suffix: 1_jd | ||
190 | args: -eps_type jd -eps_jd_blocksize 3 | ||
191 | requires: !__float128 | ||
192 | test: | ||
193 | suffix: 1_gd | ||
194 | args: -eps_type gd -eps_gd_blocksize 3 -eps_tol 1e-8 | ||
195 | test: | ||
196 | suffix: 1_gd2 | ||
197 | args: -eps_type gd -eps_gd_double_expansion | ||
198 | test: | ||
199 | suffix: 1_primme | ||
200 | args: -eps_type primme -eps_conv_abs -eps_largest_magnitude | ||
201 | requires: primme | ||
202 | |||
203 | testset: | ||
204 | args: -eps_nev 4 -eps_smallest_real -eps_max_it 600 | ||
205 | output_file: output/test8_2.out | ||
206 | test: | ||
207 | suffix: 2 | ||
208 | args: -eps_type {{rqcg lobpcg}} | ||
209 | test: | ||
210 | suffix: 2_lanczos | ||
211 | args: -eps_type lanczos -eps_lanczos_reorthog local | ||
212 | test: | ||
213 | suffix: 2_arpack | ||
214 | args: -eps_type arpack -eps_ncv 6 | ||
215 | requires: arpack !single | ||
216 | test: | ||
217 | suffix: 2_primme | ||
218 | args: -eps_type primme -eps_conv_abs -eps_primme_method lobpcg_orthobasisw -eps_ncv 24 | ||
219 | requires: primme | ||
220 | |||
221 | testset: | ||
222 | args: -eps_nev 12 -eps_mpd 9 -eps_smallest_real -eps_max_it 1000 | ||
223 | output_file: output/test8_3.out | ||
224 | test: | ||
225 | suffix: 3_rqcg | ||
226 | args: -eps_type rqcg | ||
227 | test: | ||
228 | suffix: 3_lanczos | ||
229 | args: -eps_type lanczos -eps_lanczos_reorthog local | ||
230 | test: | ||
231 | suffix: 3_lobpcg | ||
232 | args: -eps_type lobpcg -eps_lobpcg_blocksize 3 -eps_lobpcg_locking 0 -st_ksp_type preonly -st_pc_type jacobi -eps_lobpcg_checkprecond | ||
233 | requires: !__float128 | ||
234 | test: | ||
235 | suffix: 3_lobpcg_quad | ||
236 | args: -eps_type lobpcg -eps_lobpcg_blocksize 3 -eps_lobpcg_locking 0 -st_ksp_type preonly -st_pc_type jacobi -eps_tol 1e-25 | ||
237 | requires: __float128 | ||
238 | TEST*/ | ||
239 |