GCC Code Coverage Report


Directory: ./
File: src/svd/interface/svdsolve.c
Date: 2026-05-12 09:16:45
Exec Total Coverage
Lines: 258 268 96.3%
Functions: 11 11 100.0%
Branches: 729 1350 54.0%

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 SVD routines related to the solution process
12 */
13
14 #include <slepc/private/svdimpl.h> /*I "slepcsvd.h" I*/
15 #include <slepc/private/bvimpl.h>
16
17 /*
18 SVDComputeVectors_Left - Compute left singular vectors as U=A*V.
19 Only done if the leftbasis flag is false. Assumes V is available.
20 */
21 556 PetscErrorCode SVDComputeVectors_Left(SVD svd)
22 {
23 556 Vec tl,omega2,u,v,w;
24 556 PetscInt i,oldsize;
25 556 VecType vtype;
26 556 const PetscScalar* varray;
27
28
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
556 PetscFunctionBegin;
29
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
556 if (!svd->leftbasis) {
30 /* generate left singular vectors on U */
31
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.
396 if (!svd->U) PetscCall(SVDGetBV(svd,NULL,&svd->U));
32
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.
396 PetscCall(BVGetSizes(svd->U,NULL,NULL,&oldsize));
33
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
396 if (!oldsize) {
34
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.
346 if (!((PetscObject)svd->U)->type_name) PetscCall(BVSetType(svd->U,((PetscObject)svd->V)->type_name));
35
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.
346 PetscCall(MatCreateVecsEmpty(svd->A,NULL,&tl));
36
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.
346 PetscCall(BVSetSizesFromVec(svd->U,tl,svd->ncv));
37
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.
346 PetscCall(VecDestroy(&tl));
38 }
39
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.
396 PetscCall(BVSetActiveColumns(svd->V,0,svd->nconv));
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.
396 PetscCall(BVSetActiveColumns(svd->U,0,svd->nconv));
41
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.
396 if (!svd->ishyperbolic) PetscCall(BVMatMult(svd->V,svd->A,svd->U));
42
2/2
✓ Branch 0 taken 6 times.
✓ Branch 1 taken 10 times.
79 else if (svd->swapped) { /* compute right singular vectors as V=A'*Omega*U */
43
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
18 PetscCall(MatCreateVecs(svd->A,&w,NULL));
44
2/2
✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
534 for (i=0;i<svd->nconv;i++) {
45
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
516 PetscCall(BVGetColumn(svd->V,i,&v));
46
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
516 PetscCall(BVGetColumn(svd->U,i,&u));
47
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
516 PetscCall(VecPointwiseMult(w,v,svd->omega));
48
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
516 PetscCall(MatMult(svd->A,w,u));
49
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
516 PetscCall(BVRestoreColumn(svd->V,i,&v));
50
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
516 PetscCall(BVRestoreColumn(svd->U,i,&u));
51 }
52
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
18 PetscCall(VecDestroy(&w));
53 } else { /* compute left singular vectors as usual U=A*V, and set-up Omega-orthogonalization of U */
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.
61 PetscCall(BV_SetMatrixDiagonal(svd->U,svd->omega,svd->A));
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.
61 PetscCall(BVMatMult(svd->V,svd->A,svd->U));
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.
396 PetscCall(BVOrthogonalize(svd->U,NULL));
58
4/4
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 6 times.
396 if (svd->ishyperbolic && !svd->swapped) { /* store signature after Omega-orthogonalization */
59
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.
61 PetscCall(MatGetVecType(svd->A,&vtype));
60
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.
61 PetscCall(VecCreate(PETSC_COMM_SELF,&omega2));
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.
61 PetscCall(VecSetSizes(omega2,svd->nconv,svd->nconv));
62
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.
61 PetscCall(VecSetType(omega2,vtype));
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.
61 PetscCall(BVGetSignature(svd->U,omega2));
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.
61 PetscCall(VecGetArrayRead(omega2,&varray));
65
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
1869 for (i=0;i<svd->nconv;i++) {
66 1808 svd->sign[i] = PetscRealPart(varray[i]);
67
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
1808 if (PetscRealPart(varray[i])<0.0) PetscCall(BVScaleColumn(svd->U,i,-1.0));
68 }
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.
61 PetscCall(VecRestoreArrayRead(omega2,&varray));
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.
61 PetscCall(VecDestroy(&omega2));
71 }
72 }
73
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
112 PetscFunctionReturn(PETSC_SUCCESS);
74 }
75
76 24274 PetscErrorCode SVDComputeVectors(SVD svd)
77 {
78
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
24274 PetscFunctionBegin;
79
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
24274 SVDCheckSolved(svd,1);
80
8/10
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
24274 if (svd->state==SVD_STATE_SOLVED) PetscTryTypeMethod(svd,computevectors);
81 24274 svd->state = SVD_STATE_VECTORS;
82
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
24274 PetscFunctionReturn(PETSC_SUCCESS);
83 }
84
85 /*@
86 SVDSolve - Solves the singular value problem.
87
88 Collective
89
90 Input Parameter:
91 . svd - the singular value solver context
92
93 Options Database Keys:
94 + -svd_view - print information about the solver used
95 . -svd_view_mat0 - view the first matrix ($A$)
96 . -svd_view_mat1 - view the second matrix ($B$)
97 . -svd_view_signature - view the signature matrix ($\Omega$)
98 . -svd_view_vectors - view the computed singular vectors
99 . -svd_view_values - view the computed singular values
100 . -svd_converged_reason - print reason for convergence/divergence, and number of iterations
101 . -svd_error_absolute - print absolute errors of each singular triplet
102 . -svd_error_relative - print relative errors of each singular triplet
103 - -svd_error_norm - print errors relative to the matrix norms of each singular triplet
104
105 Notes:
106 The problem matrices are specified with `SVDSetOperators()`.
107
108 `SVDSolve()` will return without generating an error regardless of whether
109 all requested solutions were computed or not. Call `SVDGetConverged()` to get the
110 actual number of computed solutions, and `SVDGetConvergedReason()` to determine if
111 the solver converged or failed and why.
112
113 All the command-line options listed above admit an optional argument specifying
114 the viewer type and options. For instance, use `-svd_view_mat0 binary:amatrix.bin`
115 to save the $A$ matrix to a binary file, `-svd_view_values draw` to draw the computed
116 singular values graphically, or `-svd_error_relative :myerr.m:ascii_matlab` to save
117 the errors in a file that can be executed in Matlab.
118 See `PetscObjectViewFromOptions()` for more details.
119
120 Level: beginner
121
122 .seealso: [](ch:svd), `SVDCreate()`, `SVDSetUp()`, `SVDDestroy()`, `SVDSetOperators()`, `SVDGetConverged()`, `SVDGetConvergedReason()`
123 @*/
124 2766 PetscErrorCode SVDSolve(SVD svd)
125 {
126 2766 PetscInt i,m,n,*workperm;
127
128
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
2766 PetscFunctionBegin;
129
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
2766 PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
130
2/14
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 2 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.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
2766 if (svd->state>=SVD_STATE_SOLVED) PetscFunctionReturn(PETSC_SUCCESS);
131
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.
2766 PetscCall(PetscLogEventBegin(SVD_Solve,svd,0,0,0));
132
133 /* call setup */
134
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.
2766 PetscCall(SVDSetUp(svd));
135
136 /* safeguard for matrices with zero rows or columns */
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.
2766 PetscCall(MatGetSize(svd->OP,&m,&n));
138
2/4
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
2766 if (m == 0 || n == 0) {
139 svd->nconv = 0;
140 svd->reason = SVD_CONVERGED_TOL;
141 svd->state = SVD_STATE_SOLVED;
142 PetscFunctionReturn(PETSC_SUCCESS);
143 }
144
145 2766 svd->its = 0;
146 2766 svd->nconv = 0;
147
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
61199 for (i=0;i<svd->ncv;i++) {
148 58433 svd->sigma[i] = 0.0;
149 58433 svd->errest[i] = 0.0;
150 58433 svd->perm[i] = i;
151 }
152
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.
2766 PetscCall(SVDViewFromOptions(svd,NULL,"-svd_view_pre"));
153
154
3/4
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
2766 switch (svd->problem_type) {
155 1519 case SVD_STANDARD:
156
5/10
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
1519 PetscUseTypeMethod(svd,solve);
157 break;
158 985 case SVD_GENERALIZED:
159
5/10
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
985 PetscUseTypeMethod(svd,solveg);
160 break;
161 262 case SVD_HYPERBOLIC:
162
5/10
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
262 PetscUseTypeMethod(svd,solveh);
163 break;
164 }
165 2766 svd->state = SVD_STATE_SOLVED;
166
167 /* sort singular triplets */
168
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.
2766 if (svd->which == SVD_SMALLEST) PetscCall(PetscSortRealWithPermutation(svd->nconv,svd->sigma,svd->perm));
169 else {
170
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.
2446 PetscCall(PetscMalloc1(svd->nconv,&workperm));
171
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
35507 for (i=0;i<svd->nconv;i++) workperm[i] = i;
172
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.
2446 PetscCall(PetscSortRealWithPermutation(svd->nconv,svd->sigma,workperm));
173
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
35507 for (i=0;i<svd->nconv;i++) svd->perm[i] = workperm[svd->nconv-i-1];
174
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.
2446 PetscCall(PetscFree(workperm));
175 }
176
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.
2766 PetscCall(PetscLogEventEnd(SVD_Solve,svd,0,0,0));
177
178 /* various viewers */
179
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.
2766 PetscCall(SVDViewFromOptions(svd,NULL,"-svd_view"));
180
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.
2766 PetscCall(SVDConvergedReasonViewFromOptions(svd));
181
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.
2766 PetscCall(SVDErrorViewFromOptions(svd));
182
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.
2766 PetscCall(SVDValuesViewFromOptions(svd));
183
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.
2766 PetscCall(SVDVectorsViewFromOptions(svd));
184
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.
2766 PetscCall(MatViewFromOptions(svd->OP,(PetscObject)svd,"-svd_view_mat0"));
185
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.
2766 if (svd->isgeneralized) PetscCall(MatViewFromOptions(svd->OPb,(PetscObject)svd,"-svd_view_mat1"));
186
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.
2766 if (svd->ishyperbolic) PetscCall(VecViewFromOptions(svd->omega,(PetscObject)svd,"-svd_view_signature"));
187
188 /* Remove the initial subspaces */
189 2766 svd->nini = 0;
190 2766 svd->ninil = 0;
191
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
2766 PetscFunctionReturn(PETSC_SUCCESS);
192 }
193
194 /*@
195 SVDGetIterationNumber - Gets the current iteration number. If the
196 call to `SVDSolve()` is complete, then it returns the number of iterations
197 carried out by the solution method.
198
199 Not Collective
200
201 Input Parameter:
202 . svd - the singular value solver context
203
204 Output Parameter:
205 . its - number of iterations
206
207 Note:
208 During the $i$-th iteration this call returns $i-1$. If `SVDSolve()` is
209 complete, then parameter `its` contains either the iteration number at
210 which convergence was successfully reached, or failure was detected.
211 Call `SVDGetConvergedReason()` to determine if the solver converged or
212 failed and why.
213
214 Level: intermediate
215
216 .seealso: [](ch:svd), `SVDGetConvergedReason()`, `SVDSetTolerances()`
217 @*/
218 718 PetscErrorCode SVDGetIterationNumber(SVD svd,PetscInt *its)
219 {
220
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
718 PetscFunctionBegin;
221
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
718 PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
222
2/8
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
718 PetscAssertPointer(its,2);
223 718 *its = svd->its;
224
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
718 PetscFunctionReturn(PETSC_SUCCESS);
225 }
226
227 /*@
228 SVDGetConvergedReason - Gets the reason why the `SVDSolve()` iteration was
229 stopped.
230
231 Not Collective
232
233 Input Parameter:
234 . svd - the singular value solver context
235
236 Output Parameter:
237 . reason - negative value indicates diverged, positive value converged, see
238 `SVDConvergedReason` for the possible values
239
240 Options Database Key:
241 . -svd_converged_reason - print reason for convergence/divergence, and number of iterations
242
243 Note:
244 If this routine is called before or doing the `SVDSolve()` the value of
245 `SVD_CONVERGED_ITERATING` is returned.
246
247 Level: intermediate
248
249 .seealso: [](ch:svd), `SVDSetTolerances()`, `SVDSolve()`, `SVDConvergedReason`
250 @*/
251 171 PetscErrorCode SVDGetConvergedReason(SVD svd,SVDConvergedReason *reason)
252 {
253
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
171 PetscFunctionBegin;
254
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
171 PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
255
2/8
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
171 PetscAssertPointer(reason,2);
256
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
171 SVDCheckSolved(svd,1);
257 171 *reason = svd->reason;
258
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
171 PetscFunctionReturn(PETSC_SUCCESS);
259 }
260
261 /*@
262 SVDGetConverged - Gets the number of converged singular values.
263
264 Not Collective
265
266 Input Parameter:
267 . svd - the singular value solver context
268
269 Output Parameter:
270 . nconv - number of converged singular values
271
272 Notes:
273 This function should be called after `SVDSolve()` has finished.
274
275 The value `nconv` may be different from the number of requested solutions
276 `nsv`, but not larger than `ncv`, see `SVDSetDimensions()`.
277
278 Level: beginner
279
280 .seealso: [](ch:svd), `SVDSetDimensions()`, `SVDSolve()`, `SVDGetSingularTriplet()`
281 @*/
282 1215 PetscErrorCode SVDGetConverged(SVD svd,PetscInt *nconv)
283 {
284
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
1215 PetscFunctionBegin;
285
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
1215 PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
286
2/8
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
1215 PetscAssertPointer(nconv,2);
287
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
1215 SVDCheckSolved(svd,1);
288 1215 *nconv = svd->nconv;
289
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
1215 PetscFunctionReturn(PETSC_SUCCESS);
290 }
291
292 /*@
293 SVDGetSingularTriplet - Gets the `i`-th triplet of the singular value decomposition
294 as computed by `SVDSolve()`. The solution consists in the singular value and its left
295 and right singular vectors.
296
297 Collective
298
299 Input Parameters:
300 + svd - the singular value solver context
301 - i - index of the solution
302
303 Output Parameters:
304 + sigma - singular value
305 . u - left singular vector
306 - v - right singular vector
307
308 Note:
309 Both `u` or `v` can be `NULL` if singular vectors are not required.
310 Otherwise, the caller must provide valid `Vec` objects, i.e.,
311 they must be created by the calling program with e.g. `MatCreateVecs()`.
312
313 The index `i` should be a value between 0 and `nconv`-1 (see `SVDGetConverged()`).
314 Singular triplets are indexed according to the ordering criterion established
315 with `SVDSetWhichSingularTriplets()`.
316
317 In the case of GSVD, the solution consists in three vectors $u$, $v$, $x$ that are
318 returned as follows. Vector $x$ is returned in the right singular vector
319 (argument `v`) and has length equal to the number of columns of $A$ and $B$.
320 The other two vectors are returned stacked on top of each other $[u^*,v^*]^*$ in
321 the left singular vector argument `u`, with length equal to $m+n$ (number of rows
322 of $A$ plus number of rows of $B$).
323
324 Level: beginner
325
326 .seealso: [](ch:svd), `SVDSolve()`, `SVDGetConverged()`, `SVDSetWhichSingularTriplets()`
327 @*/
328 35257 PetscErrorCode SVDGetSingularTriplet(SVD svd,PetscInt i,PetscReal *sigma,Vec u,Vec v)
329 {
330 35257 PetscInt M,N;
331 35257 Vec w;
332
333
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
35257 PetscFunctionBegin;
334
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
35257 PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
335
27/62
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ 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 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
35257 PetscValidLogicalCollectiveInt(svd,i,2);
336
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
35257 SVDCheckSolved(svd,1);
337
17/46
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✓ Branch 22 taken 2 times.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✓ Branch 26 taken 2 times.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✓ Branch 30 taken 2 times.
✗ Branch 31 not taken.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✗ Branch 34 not taken.
✓ Branch 35 taken 2 times.
✗ Branch 36 not taken.
✗ Branch 37 not taken.
✗ Branch 38 not taken.
✓ Branch 39 taken 2 times.
✗ Branch 40 not taken.
✗ Branch 41 not taken.
✗ Branch 42 not taken.
✓ Branch 43 taken 2 times.
✗ Branch 44 not taken.
✗ Branch 45 not taken.
35257 if (u) { PetscValidHeaderSpecific(u,VEC_CLASSID,4); PetscCheckSameComm(svd,1,u,4); }
338
17/46
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✓ Branch 22 taken 2 times.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✓ Branch 26 taken 2 times.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✓ Branch 30 taken 2 times.
✗ Branch 31 not taken.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✗ Branch 34 not taken.
✓ Branch 35 taken 2 times.
✗ Branch 36 not taken.
✗ Branch 37 not taken.
✗ Branch 38 not taken.
✓ Branch 39 taken 2 times.
✗ Branch 40 not taken.
✗ Branch 41 not taken.
✗ Branch 42 not taken.
✓ Branch 43 taken 2 times.
✗ Branch 44 not taken.
✗ Branch 45 not taken.
35257 if (v) { PetscValidHeaderSpecific(v,VEC_CLASSID,5); PetscCheckSameComm(svd,1,v,5); }
339
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
35257 PetscCheck(i>=0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative");
340
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
35257 PetscCheck(i<svd->nconv,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see SVDGetConverged()");
341
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
35257 if (sigma) *sigma = svd->sigma[svd->perm[i]];
342
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
35257 if (u || v) {
343
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
24264 if (!svd->isgeneralized) {
344
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.
18217 PetscCall(MatGetSize(svd->OP,&M,&N));
345
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
18217 if (M<N) { w = u; u = v; v = w; }
346 }
347
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.
24264 PetscCall(SVDComputeVectors(svd));
348
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.
24264 if (u) PetscCall(BVCopyVec(svd->U,svd->perm[i],u));
349
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.
24264 if (v) PetscCall(BVCopyVec(svd->V,svd->perm[i],v));
350 }
351
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
7801 PetscFunctionReturn(PETSC_SUCCESS);
352 }
353
354 /*
355 SVDComputeResidualNorms_Standard - Computes the norms of the left and
356 right residuals associated with the i-th computed singular triplet.
357
358 Input Parameters:
359 sigma - singular value
360 u,v - singular vectors
361 x,y - two work vectors with the same dimensions as u,v
362 */
363 5908 static PetscErrorCode SVDComputeResidualNorms_Standard(SVD svd,PetscReal sigma,Vec u,Vec v,Vec x,Vec y,PetscReal *norm1,PetscReal *norm2)
364 {
365 5908 PetscInt M,N;
366
367
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
5908 PetscFunctionBegin;
368 /* norm1 = ||A*v-sigma*u||_2 */
369
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
5908 if (norm1) {
370
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.
5908 PetscCall(MatMult(svd->OP,v,x));
371
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.
5908 PetscCall(VecAXPY(x,-sigma,u));
372
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.
5908 PetscCall(VecNorm(x,NORM_2,norm1));
373 }
374 /* norm2 = ||A^T*u-sigma*v||_2 */
375
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
5908 if (norm2) {
376
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.
5908 PetscCall(MatGetSize(svd->OP,&M,&N));
377
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.
5908 if (M<N) PetscCall(MatMult(svd->A,u,y));
378
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.
5181 else PetscCall(MatMult(svd->AT,u,y));
379
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.
5908 PetscCall(VecAXPY(y,-sigma,v));
380
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.
5908 PetscCall(VecNorm(y,NORM_2,norm2));
381 }
382
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
927 PetscFunctionReturn(PETSC_SUCCESS);
383 }
384
385 /*
386 SVDComputeResidualNorms_Generalized - In GSVD, compute the residual norms
387 norm1 = ||s^2*A'*u-c*B'*B*x||_2 and norm2 = ||c^2*B'*v-s*A'*A*x||_2.
388
389 Input Parameters:
390 sigma - singular value
391 uv - left singular vectors [u;v]
392 x - right singular vector
393 y,z - two work vectors with the same dimension as x
394 */
395 4077 static PetscErrorCode SVDComputeResidualNorms_Generalized(SVD svd,PetscReal sigma,Vec uv,Vec x,Vec y,Vec z,PetscReal *norm1,PetscReal *norm2)
396 {
397 4077 Vec u,v,ax,bx,nest,aux[2];
398 4077 PetscReal c,s;
399
400
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
4077 PetscFunctionBegin;
401
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.
4077 PetscCall(MatCreateVecs(svd->OP,NULL,&u));
402
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.
4077 PetscCall(MatCreateVecs(svd->OPb,NULL,&v));
403 4077 aux[0] = u;
404 4077 aux[1] = v;
405
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.
4077 PetscCall(VecCreateNest(PetscObjectComm((PetscObject)svd),2,NULL,aux,&nest));
406
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.
4077 PetscCall(VecCopy(uv,nest));
407
408 4077 s = 1.0/PetscSqrtReal(1.0+sigma*sigma);
409 4077 c = sigma*s;
410
411 /* norm1 = ||s^2*A'*u-c*B'*B*x||_2 */
412
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
4077 if (norm1) {
413
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.
4077 PetscCall(VecDuplicate(v,&bx));
414
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.
4077 PetscCall(MatMultHermitianTranspose(svd->OP,u,z));
415
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.
4077 PetscCall(MatMult(svd->OPb,x,bx));
416
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.
4077 PetscCall(MatMultHermitianTranspose(svd->OPb,bx,y));
417
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.
4077 PetscCall(VecAXPBY(y,s*s,-c,z));
418
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.
4077 PetscCall(VecNorm(y,NORM_2,norm1));
419
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.
4077 PetscCall(VecDestroy(&bx));
420 }
421 /* norm2 = ||c^2*B'*v-s*A'*A*x||_2 */
422
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
4077 if (norm2) {
423
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.
4077 PetscCall(VecDuplicate(u,&ax));
424
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.
4077 PetscCall(MatMultHermitianTranspose(svd->OPb,v,z));
425
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.
4077 PetscCall(MatMult(svd->OP,x,ax));
426
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.
4077 PetscCall(MatMultHermitianTranspose(svd->OP,ax,y));
427
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.
4077 PetscCall(VecAXPBY(y,c*c,-s,z));
428
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.
4077 PetscCall(VecNorm(y,NORM_2,norm2));
429
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.
4077 PetscCall(VecDestroy(&ax));
430 }
431
432
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.
4077 PetscCall(VecDestroy(&v));
433
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.
4077 PetscCall(VecDestroy(&u));
434
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.
4077 PetscCall(VecDestroy(&nest));
435
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
807 PetscFunctionReturn(PETSC_SUCCESS);
436 }
437
438 /*
439 SVDComputeResidualNorms_Hyperbolic - Computes the norms of the left and
440 right residuals associated with the i-th computed singular triplet.
441
442 Input Parameters:
443 sigma - singular value
444 sign - corresponding element of the signature Omega2
445 u,v - singular vectors
446 x,y,z - three work vectors with the same dimensions as u,v,u
447 */
448 4166 static PetscErrorCode SVDComputeResidualNorms_Hyperbolic(SVD svd,PetscReal sigma,PetscReal sign,Vec u,Vec v,Vec x,Vec y,Vec z,PetscReal *norm1,PetscReal *norm2)
449 {
450 4166 PetscInt M,N;
451
452
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
4166 PetscFunctionBegin;
453 /* norm1 = ||A*v-sigma*u||_2 */
454
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
4166 if (norm1) {
455
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.
4166 PetscCall(MatMult(svd->OP,v,x));
456
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.
4166 PetscCall(VecAXPY(x,-sigma,u));
457
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.
4166 PetscCall(VecNorm(x,NORM_2,norm1));
458 }
459 /* norm2 = ||A^T*Omega*u-sigma*sign*v||_2 */
460
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
4166 if (norm2) {
461
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.
4166 PetscCall(MatGetSize(svd->OP,&M,&N));
462
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.
4166 PetscCall(VecPointwiseMult(z,u,svd->omega));
463
6/8
✓ Branch 0 taken 6 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
4166 if (M<N) PetscCall(MatMult(svd->A,z,y));
464
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.
3878 else PetscCall(MatMult(svd->AT,z,y));
465
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.
4166 PetscCall(VecAXPY(y,-sigma*sign,v));
466
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.
4166 PetscCall(VecNorm(y,NORM_2,norm2));
467 }
468
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
1274 PetscFunctionReturn(PETSC_SUCCESS);
469 }
470
471 /*@
472 SVDComputeError - Computes the error (based on the residual norm) associated
473 with the `i`-th singular triplet.
474
475 Collective
476
477 Input Parameters:
478 + svd - the singular value solver context
479 . i - the solution index
480 - type - the type of error to compute, see `SVDErrorType`
481
482 Output Parameter:
483 . error - the error
484
485 Notes:
486 The error can be computed in various ways, all of them based on the residual
487 norm obtained as $\sqrt{\eta_1^2+\eta_2^2}$ with $\eta_1 = \|Av-\sigma u\|_2$ and
488 $\eta_2 = \|A^*u-\sigma v\|_2$, where $(\sigma,u,v)$ is the approximate singular
489 triplet.
490
491 In the case of the GSVD, the two components of the residual norm are
492 $\eta_1 = \|s^2 A^*u-cB^*Bx\|_2$ and $\eta_2 = ||c^2 B^*v-sA^*Ax||_2$, where
493 $(\sigma,u,v,x)$ is the approximate generalized singular quadruple, with
494 $\sigma=c/s$.
495
496 Level: beginner
497
498 .seealso: [](ch:svd), `SVDErrorType`, `SVDSolve()`
499 @*/
500 14151 PetscErrorCode SVDComputeError(SVD svd,PetscInt i,SVDErrorType type,PetscReal *error)
501 {
502 14151 PetscReal sigma,norm1,norm2,c,s;
503 14151 Vec u=NULL,v=NULL,x=NULL,y=NULL,z=NULL;
504 14151 PetscReal vecnorm=1.0;
505
506
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
14151 PetscFunctionBegin;
507
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
14151 PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
508
27/62
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ 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 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
14151 PetscValidLogicalCollectiveInt(svd,i,2);
509
27/62
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ 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 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
14151 PetscValidLogicalCollectiveEnum(svd,type,3);
510
2/8
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
14151 PetscAssertPointer(error,4);
511
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
14151 SVDCheckSolved(svd,1);
512
513 /* allocate work vectors */
514
3/4
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
14151 switch (svd->problem_type) {
515 5908 case SVD_STANDARD:
516
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.
5908 PetscCall(SVDSetWorkVecs(svd,2,2));
517 5908 u = svd->workl[0];
518 5908 v = svd->workr[0];
519 5908 x = svd->workl[1];
520 5908 y = svd->workr[1];
521 5908 break;
522 4077 case SVD_GENERALIZED:
523
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.
4077 PetscCall(SVDSetWorkVecs(svd,1,3));
524 4077 u = svd->workl[0];
525 4077 v = svd->workr[0];
526 4077 x = svd->workr[1];
527 4077 y = svd->workr[2];
528 4077 break;
529 4166 case SVD_HYPERBOLIC:
530
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.
4166 PetscCall(SVDSetWorkVecs(svd,3,2));
531 4166 u = svd->workl[0];
532 4166 v = svd->workr[0];
533 4166 x = svd->workl[1];
534 4166 y = svd->workr[1];
535 4166 z = svd->workl[2];
536 4166 break;
537 }
538
539 /* compute residual norm */
540
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.
14151 PetscCall(SVDGetSingularTriplet(svd,i,&sigma,u,v));
541
3/4
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
14151 switch (svd->problem_type) {
542 5908 case SVD_STANDARD:
543
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.
5908 PetscCall(SVDComputeResidualNorms_Standard(svd,sigma,u,v,x,y,&norm1,&norm2));
544 break;
545 4077 case SVD_GENERALIZED:
546
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.
4077 PetscCall(SVDComputeResidualNorms_Generalized(svd,sigma,u,v,x,y,&norm1,&norm2));
547 break;
548 4166 case SVD_HYPERBOLIC:
549
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.
4166 PetscCall(SVDComputeResidualNorms_Hyperbolic(svd,sigma,svd->sign[svd->perm[i]],u,v,x,y,z,&norm1,&norm2));
550 break;
551 }
552 14151 *error = SlepcAbs(norm1,norm2);
553
554 /* compute 2-norm of eigenvector of the cyclic form */
555
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
14151 if (type!=SVD_ERROR_ABSOLUTE) {
556
3/4
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
14017 switch (svd->problem_type) {
557 5774 case SVD_STANDARD:
558 5774 vecnorm = PETSC_SQRT2;
559 5774 break;
560 4077 case SVD_GENERALIZED:
561
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.
4077 PetscCall(VecNorm(v,NORM_2,&vecnorm));
562 4077 vecnorm = PetscSqrtReal(1.0+vecnorm*vecnorm);
563 4077 break;
564 4166 case SVD_HYPERBOLIC:
565
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.
4166 PetscCall(VecNorm(u,NORM_2,&vecnorm));
566 4166 vecnorm = PetscSqrtReal(1.0+vecnorm*vecnorm);
567 4166 break;
568 }
569 13 }
570
571 /* compute error */
572
3/4
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
14030 switch (type) {
573 case SVD_ERROR_ABSOLUTE:
574 break;
575 5700 case SVD_ERROR_RELATIVE:
576
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
5700 if (svd->isgeneralized) {
577 s = 1.0/PetscSqrtReal(1.0+sigma*sigma);
578 c = sigma*s;
579 norm1 /= c*vecnorm;
580 norm2 /= s*vecnorm;
581 *error = PetscMax(norm1,norm2);
582 5700 } else *error /= sigma*vecnorm;
583 break;
584 8317 case SVD_ERROR_NORM:
585
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.
8317 if (!svd->nrma) PetscCall(MatNorm(svd->OP,NORM_INFINITY,&svd->nrma));
586
8/10
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
8317 if (svd->isgeneralized && !svd->nrmb) PetscCall(MatNorm(svd->OPb,NORM_INFINITY,&svd->nrmb));
587
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
8317 *error /= PetscMax(svd->nrma,svd->nrmb)*vecnorm;
588 8317 break;
589 default:
590 SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
591 }
592
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
3008 PetscFunctionReturn(PETSC_SUCCESS);
593 }
594