GCC Code Coverage Report


Directory: ./
File: src/svd/tests/test3.c
Date: 2025-10-04 04:19:13
Exec Total Coverage
Lines: 50 51 98.0%
Functions: 1 1 100.0%
Branches: 154 236 65.3%

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[] = "Test SVD with user-provided initial vectors.\n\n"
12 "The command line options are:\n"
13 " -n <n>, where <n> = row dimension.\n"
14 " -m <m>, where <m> = column dimension.\n\n";
15
16 #include <slepcsvd.h>
17
18 /*
19 This example computes the singular values of a rectangular nxm Grcar matrix:
20
21 | 1 1 1 1 |
22 | -1 1 1 1 1 |
23 | -1 1 1 1 1 |
24 A = | . . . . . |
25 | . . . . . |
26 | -1 1 1 1 1 |
27 | -1 1 1 1 |
28
29 */
30
31 309 int main(int argc,char **argv)
32 {
33 309 Mat A; /* Grcar matrix */
34 309 SVD svd; /* singular value solver context */
35 309 Vec v0,w0; /* initial vectors */
36 309 Vec *U,*V;
37 309 PetscInt N=35,M=30,Istart,Iend,i,col[5],nconv;
38 309 PetscScalar value[] = { -1, 1, 1, 1, 1 };
39 309 PetscReal lev1=0.0,lev2=0.0,tol=PETSC_SMALL;
40 309 PetscBool skiporth=PETSC_FALSE;
41
42
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
309 PetscFunctionBeginUser;
43
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.
309 PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
44
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.
309 PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
45
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.
309 PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&M,NULL));
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.
309 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSVD of a rectangular Grcar matrix, %" PetscInt_FMT "x%" PetscInt_FMT "\n\n",N,M));
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.
309 PetscCall(PetscOptionsGetBool(NULL,NULL,"-skiporth",&skiporth,NULL));
48
49 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50 Generate the matrix
51 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
52
53
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.
309 PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
54
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.
309 PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,M));
55
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.
309 PetscCall(MatSetFromOptions(A));
56
57
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.
309 PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
58
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
10774 for (i=Istart;i<Iend;i++) {
59 10465 col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
60
6/8
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
10465 if (i==0) PetscCall(MatSetValues(A,1,&i,PetscMin(4,M-i+1),col+1,value+1,INSERT_VALUES));
61
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.
10465 else PetscCall(MatSetValues(A,1,&i,PetscMin(5,M-i+1),col,value,INSERT_VALUES));
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.
309 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
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.
309 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
65
66 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67 Create the SVD context and solve the problem
68 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
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.
309 PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
71
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.
309 PetscCall(SVDSetOperators(svd,A,NULL));
72
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.
309 PetscCall(SVDSetFromOptions(svd));
73
74 /*
75 Set the initial vectors. This is optional, if not done the initial
76 vectors are set to random values
77 */
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.
309 PetscCall(MatCreateVecs(A,&v0,&w0));
79
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.
309 PetscCall(VecSet(v0,1.0));
80
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.
309 PetscCall(VecSet(w0,1.0));
81
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.
309 PetscCall(SVDSetInitialSpaces(svd,1,&v0,1,&w0));
82
83 /*
84 Compute solution
85 */
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.
309 PetscCall(SVDSolve(svd));
87
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.
309 PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL));
88
89 /*
90 Check orthonormality of computed singular vectors
91 */
92
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.
309 PetscCall(SVDGetConverged(svd,&nconv));
93
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
309 if (nconv>1) {
94
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.
309 PetscCall(VecDuplicateVecs(w0,nconv,&U));
95
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.
309 PetscCall(VecDuplicateVecs(v0,nconv,&V));
96
7/8
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
3172 for (i=0;i<nconv;i++) PetscCall(SVDGetSingularTriplet(svd,i,NULL,U[i],V[i]));
97
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
309 if (!skiporth) {
98
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.
309 PetscCall(VecCheckOrthonormality(U,nconv,NULL,nconv,NULL,NULL,&lev1));
99
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.
309 PetscCall(VecCheckOrthonormality(V,nconv,NULL,nconv,NULL,NULL,&lev2));
100 }
101
5/8
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
309 if (lev1+lev2<20*tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality below the tolerance\n"));
102 else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g (U) %g (V)\n",(double)lev1,(double)lev2));
103
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.
309 PetscCall(VecDestroyVecs(nconv,&U));
104
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.
309 PetscCall(VecDestroyVecs(nconv,&V));
105 }
106
107 /*
108 Free work space
109 */
110
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.
309 PetscCall(VecDestroy(&v0));
111
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.
309 PetscCall(VecDestroy(&w0));
112
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.
309 PetscCall(SVDDestroy(&svd));
113
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.
309 PetscCall(MatDestroy(&A));
114
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.
309 PetscCall(SlepcFinalize());
115 return 0;
116 }
117
118 /*TEST
119
120 testset:
121 args: -svd_nsv 4
122 output_file: output/test3_1.out
123 filter: sed -e "s/22176/22175/" | sed -e "s/21798/21797/" | sed -e "s/16826/16825/" | sed -e "s/15129/15128/" | sed -e "s/22200/22201/"
124 test:
125 suffix: 1_lanczos
126 args: -svd_type lanczos
127 test:
128 suffix: 1_lanczos_one
129 args: -svd_type lanczos -svd_lanczos_oneside
130 test:
131 suffix: 1_trlanczos
132 args: -svd_type trlanczos -svd_trlanczos_locking {{0 1}}
133 test:
134 suffix: 1_trlanczos_one
135 args: -svd_type trlanczos -svd_trlanczos_oneside
136 test:
137 suffix: 1_trlanczos_one_mgs
138 args: -svd_type trlanczos -svd_trlanczos_oneside -bv_orthog_type mgs
139 test:
140 suffix: 1_trlanczos_one_always
141 args: -svd_type trlanczos -svd_trlanczos_oneside -bv_orthog_refine always
142 test:
143 suffix: 1_cross
144 args: -svd_type cross
145 test:
146 suffix: 1_cross_exp
147 args: -svd_type cross -svd_cross_explicitmatrix
148 test:
149 suffix: 1_cyclic
150 args: -svd_type cyclic
151 requires: !__float128
152 test:
153 suffix: 1_cyclic_exp
154 args: -svd_type cyclic -svd_cyclic_explicitmatrix
155 requires: !__float128
156 test:
157 suffix: 1_lapack
158 args: -svd_type lapack
159 test:
160 suffix: 1_randomized
161 args: -svd_type randomized
162 test:
163 suffix: 1_primme
164 args: -svd_type primme
165 requires: primme
166
167 testset:
168 args: -svd_implicittranspose -svd_nsv 4 -svd_tol 1e-5
169 output_file: output/test3_1.out
170 filter: sed -e "s/22176/22175/" | sed -e "s/21798/21797/" | sed -e "s/16826/16825/" | sed -e "s/15129/15128/" | sed -e "s/22200/22201/"
171 test:
172 suffix: 2_lanczos
173 args: -svd_type lanczos -svd_conv_norm
174 test:
175 suffix: 2_lanczos_one
176 args: -svd_type lanczos -svd_lanczos_oneside
177 test:
178 suffix: 2_trlanczos
179 args: -svd_type trlanczos
180 test:
181 suffix: 2_trlanczos_one
182 args: -svd_type trlanczos -svd_trlanczos_oneside
183 test:
184 suffix: 2_trlanczos_one_mgs
185 args: -svd_type trlanczos -svd_trlanczos_oneside -bv_orthog_type mgs
186 test:
187 suffix: 2_trlanczos_one_always
188 args: -svd_type trlanczos -svd_trlanczos_oneside -bv_orthog_refine always
189 test:
190 suffix: 2_cross
191 args: -svd_type cross
192 test:
193 suffix: 2_cross_exp
194 args: -svd_type cross -svd_cross_explicitmatrix
195 requires: !complex
196 test:
197 suffix: 2_cyclic
198 args: -svd_type cyclic -svd_tol 1e-9
199 requires: double
200 test:
201 suffix: 2_lapack
202 args: -svd_type lapack
203 test:
204 suffix: 2_randomized
205 args: -svd_type randomized
206
207 testset:
208 args: -svd_nsv 4 -mat_type aijcusparse
209 requires: cuda
210 output_file: output/test3_1.out
211 filter: sed -e "s/22176/22175/" | sed -e "s/21798/21797/" | sed -e "s/16826/16825/" | sed -e "s/15129/15128/"
212 test:
213 suffix: 3_cuda_lanczos
214 args: -svd_type lanczos
215 test:
216 suffix: 3_cuda_lanczos_one
217 args: -svd_type lanczos -svd_lanczos_oneside
218 test:
219 suffix: 3_cuda_trlanczos
220 args: -svd_type trlanczos
221 test:
222 suffix: 3_cuda_trlanczos_one
223 args: -svd_type trlanczos -svd_trlanczos_oneside
224 test:
225 suffix: 3_cuda_trlanczos_one_mgs
226 args: -svd_type trlanczos -svd_trlanczos_oneside -bv_orthog_type mgs
227 test:
228 suffix: 3_cuda_trlanczos_one_always
229 args: -svd_type trlanczos -svd_trlanczos_oneside -bv_orthog_refine always
230 test:
231 suffix: 3_cuda_cross
232 args: -svd_type cross
233 test:
234 suffix: 3_cuda_cyclic
235 args: -svd_type cyclic
236 test:
237 suffix: 3_cuda_cyclic_exp
238 args: -svd_type cyclic -svd_cyclic_explicitmatrix
239 test:
240 suffix: 3_cuda_randomized
241 args: -svd_type randomized
242
243 test:
244 suffix: 4
245 args: -svd_type lapack -svd_nsv 4
246 output_file: output/test3_1.out
247 filter: sed -e "s/15129/15128/" | sed -e "s/21798/21797/"
248 nsize: 2
249
250 test:
251 suffix: 5
252 args: -svd_nsv 4 -svd_view_values draw -svd_monitor draw::draw_lg
253 requires: x
254 output_file: output/test3_1.out
255
256 testset:
257 args: -svd_nsv 4 -mat_type aijhipsparse
258 requires: hip
259 output_file: output/test3_1.out
260 filter: sed -e "s/22176/22175/" | sed -e "s/21798/21797/" | sed -e "s/16826/16825/" | sed -e "s/15129/15128/"
261 test:
262 suffix: 6_hip_lanczos
263 args: -svd_type lanczos
264 test:
265 suffix: 6_hip_lanczos_one
266 args: -svd_type lanczos -svd_lanczos_oneside
267 test:
268 suffix: 6_hip_trlanczos
269 args: -svd_type trlanczos
270 test:
271 suffix: 6_hip_trlanczos_one
272 args: -svd_type trlanczos -svd_trlanczos_oneside
273 test:
274 suffix: 6_hip_trlanczos_one_mgs
275 args: -svd_type trlanczos -svd_trlanczos_oneside -bv_orthog_type mgs
276 test:
277 suffix: 6_hip_trlanczos_one_always
278 args: -svd_type trlanczos -svd_trlanczos_oneside -bv_orthog_refine always
279 test:
280 suffix: 6_hip_cross
281 args: -svd_type cross
282 test:
283 suffix: 6_hip_cyclic
284 args: -svd_type cyclic
285 test:
286 suffix: 6_hip_cyclic_exp
287 args: -svd_type cyclic -svd_cyclic_explicitmatrix
288 test:
289 suffix: 6_hip_randomized
290 args: -svd_type randomized
291
292 TEST*/
293