GCC Code Coverage Report


Directory: ./
File: src/svd/tests/test12.c
Date: 2025-10-04 04:19:13
Exec Total Coverage
Lines: 57 59 96.6%
Functions: 2 2 100.0%
Branches: 143 250 57.2%

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[] = "SVD problem with user-defined stopping test.\n\n"
12 "The command line options are:\n"
13 " -m <m>, where <m> = matrix rows.\n"
14 " -n <n>, where <n> = matrix columns (defaults to m+2).\n\n";
15
16 #include <slepcsvd.h>
17 #include <petsctime.h>
18
19 /*
20 This example computes the singular values of a rectangular bidiagonal matrix
21
22 | 1 2 |
23 | 1 2 |
24 | 1 2 |
25 A = | . . |
26 | . . |
27 | 1 2 |
28 | 1 2 |
29 */
30
31 /*
32 Function for user-defined stopping test.
33
34 Checks that the computing time has not exceeded the deadline.
35 */
36 5204 PetscErrorCode MyStoppingTest(SVD svd,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,SVDConvergedReason *reason,void *ctx)
37 {
38 5204 PetscLogDouble now,deadline = *(PetscLogDouble*)ctx;
39
40
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
5204 PetscFunctionBeginUser;
41 /* check if usual termination conditions are met */
42
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.
5204 PetscCall(SVDStoppingBasic(svd,its,max_it,nconv,nev,reason,NULL));
43
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
5204 if (*reason==SVD_CONVERGED_ITERATING) {
44 /* check if deadline has expired */
45
2/4
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
5204 PetscCall(PetscTime(&now));
46
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
5204 if (now>deadline) *reason = SVD_CONVERGED_USER;
47 }
48
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.
369 PetscFunctionReturn(PETSC_SUCCESS);
49 }
50
51 10 int main(int argc,char **argv)
52 {
53 10 Mat A;
54 10 SVD svd;
55 10 SVDConvergedReason reason;
56 10 PetscInt m=20,n,Istart,Iend,i,col[2],nconv;
57 10 PetscReal seconds=2.5; /* maximum time allowed for computation */
58 10 PetscLogDouble deadline; /* time to abort computation */
59 10 PetscScalar value[] = { 1, 2 };
60 10 PetscBool terse,flg;
61 10 PetscViewer viewer;
62
63
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
10 PetscFunctionBeginUser;
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.
10 PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
65
66
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(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
67
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(PetscOptionsGetInt(NULL,NULL,"-n",&n,&flg));
68
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
10 if (!flg) n=m+2;
69
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(PetscPrintf(PETSC_COMM_WORLD,"\nRectangular bidiagonal matrix, m=%" PetscInt_FMT " n=%" PetscInt_FMT "\n\n",m,n));
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.
10 PetscCall(PetscOptionsGetReal(NULL,NULL,"-seconds",&seconds,NULL));
71 10 deadline = seconds;
72
2/4
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
10 PetscCall(PetscTimeAdd(&deadline));
73
74 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75 Generate the matrix
76 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
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.
10 PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
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.
10 PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n));
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.
10 PetscCall(MatSetFromOptions(A));
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.
10 PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
82
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
7510 for (i=Istart;i<Iend;i++) {
83 7500 col[0]=i; col[1]=i+1;
84
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.
7500 if (i<n-1) PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES));
85
0/8
✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
7500 else if (i==n-1) PetscCall(MatSetValue(A,i,col[0],value[0],INSERT_VALUES));
86 }
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.
10 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
88
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(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
89
90 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91 Compute singular values
92 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93
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.
10 PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
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.
10 PetscCall(SVDSetOperators(svd,A,NULL));
96
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(SVDSetWhichSingularTriplets(svd,SVD_SMALLEST));
97
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(SVDSetTolerances(svd,PETSC_CURRENT,1000));
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.
10 PetscCall(SVDSetType(svd,SVDTRLANCZOS));
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.
10 PetscCall(SVDSetStoppingTestFunction(svd,MyStoppingTest,&deadline,NULL));
100
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(SVDSetFromOptions(svd));
101
102 /* call the solver */
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.
10 PetscCall(SVDSolve(svd));
104
105 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106 Display solution and clean up
107 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108
109 /* show detailed info unless -terse option is given by user */
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.
10 PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
111
1/8
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
10 if (terse) PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
112 else {
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.
10 PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
114
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(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
115
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(SVDGetConvergedReason(svd,&reason));
116
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
10 if (reason!=SVD_CONVERGED_USER) {
117 PetscCall(SVDConvergedReasonView(svd,viewer));
118 PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,viewer));
119 } else {
120
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(SVDGetConverged(svd,&nconv));
121
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(PetscViewerASCIIPrintf(viewer,"SVD solve finished with %" PetscInt_FMT " converged eigenpairs; reason=%s\n",nconv,SVDConvergedReasons[reason]));
122 }
123
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(PetscViewerPopFormat(viewer));
124 }
125
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(SVDDestroy(&svd));
126
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(MatDestroy(&A));
127
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.
10 PetscCall(SlepcFinalize());
128 return 0;
129 }
130
131 /*TEST
132
133 test:
134 suffix: 1
135 args: -m 750 -seconds 0.1 -svd_max_it 10000
136 requires: !single
137
138 TEST*/
139