GCC Code Coverage Report


Directory: ./
File: src/svd/interface/svdsolve.c
Date: 2025-11-19 04:19:03
Exec Total Coverage
Lines: 258 269 95.9%
Functions: 11 11 100.0%
Branches: 729 1352 53.9%

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 563 PetscErrorCode SVDComputeVectors_Left(SVD svd)
22 {
23 563 Vec tl,omega2,u,v,w;
24 563 PetscInt i,oldsize;
25 563 VecType vtype;
26 563 const PetscScalar* varray;
27
28
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
563 PetscFunctionBegin;
29
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
563 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.
392 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.
392 PetscCall(BVGetSizes(svd->U,NULL,NULL,&oldsize));
33
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
392 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.
342 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.
342 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.
342 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.
342 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.
392 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.
392 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.
392 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.
392 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.
392 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.
1877 for (i=0;i<svd->nconv;i++) {
66 1816 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.
1816 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.
113 PetscFunctionReturn(PETSC_SUCCESS);
74 }
75
76 24287 PetscErrorCode SVDComputeVectors(SVD svd)
77 {
78
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
24287 PetscFunctionBegin;
79
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
24287 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.
24287 if (svd->state==SVD_STATE_SOLVED) PetscTryTypeMethod(svd,computevectors);
81 24287 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.
24287 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
119 Level: beginner
120
121 .seealso: [](ch:svd), `SVDCreate()`, `SVDSetUp()`, `SVDDestroy()`, `SVDSetOperators()`, `SVDGetConverged()`, `SVDGetConvergedReason()`
122 @*/
123 2773 PetscErrorCode SVDSolve(SVD svd)
124 {
125 2773 PetscInt i,m,n,*workperm;
126
127
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
2773 PetscFunctionBegin;
128
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.
2773 PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
129
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.
2773 if (svd->state>=SVD_STATE_SOLVED) PetscFunctionReturn(PETSC_SUCCESS);
130
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.
2773 PetscCall(PetscLogEventBegin(SVD_Solve,svd,0,0,0));
131
132 /* call setup */
133
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.
2773 PetscCall(SVDSetUp(svd));
134
135 /* safeguard for matrices with zero rows or columns */
136
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.
2773 PetscCall(MatGetSize(svd->OP,&m,&n));
137
2/4
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
2773 if (m == 0 || n == 0) {
138 svd->nconv = 0;
139 svd->reason = SVD_CONVERGED_TOL;
140 svd->state = SVD_STATE_SOLVED;
141 PetscFunctionReturn(PETSC_SUCCESS);
142 }
143
144 2773 svd->its = 0;
145 2773 svd->nconv = 0;
146
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
61325 for (i=0;i<svd->ncv;i++) {
147 58552 svd->sigma[i] = 0.0;
148 58552 svd->errest[i] = 0.0;
149 58552 svd->perm[i] = i;
150 }
151
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.
2773 PetscCall(SVDViewFromOptions(svd,NULL,"-svd_view_pre"));
152
153
3/4
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
2773 switch (svd->problem_type) {
154 1526 case SVD_STANDARD:
155
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.
1526 PetscUseTypeMethod(svd,solve);
156 break;
157 985 case SVD_GENERALIZED:
158
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);
159 break;
160 262 case SVD_HYPERBOLIC:
161
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);
162 break;
163 }
164 2773 svd->state = SVD_STATE_SOLVED;
165
166 /* sort singular triplets */
167
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.
2773 if (svd->which == SVD_SMALLEST) PetscCall(PetscSortRealWithPermutation(svd->nconv,svd->sigma,svd->perm));
168 else {
169
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.
2453 PetscCall(PetscMalloc1(svd->nconv,&workperm));
170
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
35568 for (i=0;i<svd->nconv;i++) workperm[i] = i;
171
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.
2453 PetscCall(PetscSortRealWithPermutation(svd->nconv,svd->sigma,workperm));
172
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
35568 for (i=0;i<svd->nconv;i++) svd->perm[i] = workperm[svd->nconv-i-1];
173
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.
2453 PetscCall(PetscFree(workperm));
174 }
175
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.
2773 PetscCall(PetscLogEventEnd(SVD_Solve,svd,0,0,0));
176
177 /* various viewers */
178
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.
2773 PetscCall(SVDViewFromOptions(svd,NULL,"-svd_view"));
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.
2773 PetscCall(SVDConvergedReasonViewFromOptions(svd));
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.
2773 PetscCall(SVDErrorViewFromOptions(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.
2773 PetscCall(SVDValuesViewFromOptions(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.
2773 PetscCall(SVDVectorsViewFromOptions(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.
2773 PetscCall(MatViewFromOptions(svd->OP,(PetscObject)svd,"-svd_view_mat0"));
184
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.
2773 if (svd->isgeneralized) PetscCall(MatViewFromOptions(svd->OPb,(PetscObject)svd,"-svd_view_mat1"));
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.
2773 if (svd->ishyperbolic) PetscCall(VecViewFromOptions(svd->omega,(PetscObject)svd,"-svd_view_signature"));
186
187 /* Remove the initial subspaces */
188 2773 svd->nini = 0;
189 2773 svd->ninil = 0;
190
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.
2773 PetscFunctionReturn(PETSC_SUCCESS);
191 }
192
193 /*@
194 SVDGetIterationNumber - Gets the current iteration number. If the
195 call to `SVDSolve()` is complete, then it returns the number of iterations
196 carried out by the solution method.
197
198 Not Collective
199
200 Input Parameter:
201 . svd - the singular value solver context
202
203 Output Parameter:
204 . its - number of iterations
205
206 Note:
207 During the $i$-th iteration this call returns $i-1$. If `SVDSolve()` is
208 complete, then parameter `its` contains either the iteration number at
209 which convergence was successfully reached, or failure was detected.
210 Call `SVDGetConvergedReason()` to determine if the solver converged or
211 failed and why.
212
213 Level: intermediate
214
215 .seealso: [](ch:svd), `SVDGetConvergedReason()`, `SVDSetTolerances()`
216 @*/
217 729 PetscErrorCode SVDGetIterationNumber(SVD svd,PetscInt *its)
218 {
219
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
729 PetscFunctionBegin;
220
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.
729 PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
221
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.
729 PetscAssertPointer(its,2);
222 729 *its = svd->its;
223
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.
729 PetscFunctionReturn(PETSC_SUCCESS);
224 }
225
226 /*@
227 SVDGetConvergedReason - Gets the reason why the `SVDSolve()` iteration was
228 stopped.
229
230 Not Collective
231
232 Input Parameter:
233 . svd - the singular value solver context
234
235 Output Parameter:
236 . reason - negative value indicates diverged, positive value converged, see
237 `SVDConvergedReason` for the possible values
238
239 Options Database Key:
240 . -svd_converged_reason - print reason for convergence/divergence, and number of iterations
241
242 Note:
243 If this routine is called before or doing the `SVDSolve()` the value of
244 `SVD_CONVERGED_ITERATING` is returned.
245
246 Level: intermediate
247
248 .seealso: [](ch:svd), `SVDSetTolerances()`, `SVDSolve()`, `SVDConvergedReason`
249 @*/
250 171 PetscErrorCode SVDGetConvergedReason(SVD svd,SVDConvergedReason *reason)
251 {
252
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
171 PetscFunctionBegin;
253
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);
254
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);
255
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
171 SVDCheckSolved(svd,1);
256 171 *reason = svd->reason;
257
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);
258 }
259
260 /*@
261 SVDGetConverged - Gets the number of converged singular values.
262
263 Not Collective
264
265 Input Parameter:
266 . svd - the singular value solver context
267
268 Output Parameter:
269 . nconv - number of converged singular values
270
271 Notes:
272 This function should be called after `SVDSolve()` has finished.
273
274 The value `nconv` may be different from the number of requested solutions
275 `nsv`, but not larger than `ncv`, see `SVDSetDimensions()`.
276
277 Level: beginner
278
279 .seealso: [](ch:svd), `SVDSetDimensions()`, `SVDSolve()`, `SVDGetSingularTriplet()`
280 @*/
281 1213 PetscErrorCode SVDGetConverged(SVD svd,PetscInt *nconv)
282 {
283
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
1213 PetscFunctionBegin;
284
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.
1213 PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
285
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.
1213 PetscAssertPointer(nconv,2);
286
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
1213 SVDCheckSolved(svd,1);
287 1213 *nconv = svd->nconv;
288
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.
1213 PetscFunctionReturn(PETSC_SUCCESS);
289 }
290
291 /*@
292 SVDGetSingularTriplet - Gets the `i`-th triplet of the singular value decomposition
293 as computed by `SVDSolve()`. The solution consists in the singular value and its left
294 and right singular vectors.
295
296 Collective
297
298 Input Parameters:
299 + svd - the singular value solver context
300 - i - index of the solution
301
302 Output Parameters:
303 + sigma - singular value
304 . u - left singular vector
305 - v - right singular vector
306
307 Note:
308 Both `u` or `v` can be `NULL` if singular vectors are not required.
309 Otherwise, the caller must provide valid `Vec` objects, i.e.,
310 they must be created by the calling program with e.g. `MatCreateVecs()`.
311
312 The index `i` should be a value between 0 and `nconv`-1 (see `SVDGetConverged()`).
313 Singular triplets are indexed according to the ordering criterion established
314 with `SVDSetWhichSingularTriplets()`.
315
316 In the case of GSVD, the solution consists in three vectors $u$, $v$, $x$ that are
317 returned as follows. Vector $x$ is returned in the right singular vector
318 (argument `v`) and has length equal to the number of columns of $A$ and $B$.
319 The other two vectors are returned stacked on top of each other $[u^*,v^*]^*$ in
320 the left singular vector argument `u`, with length equal to $m+n$ (number of rows
321 of $A$ plus number of rows of $B$).
322
323 Level: beginner
324
325 .seealso: [](ch:svd), `SVDSolve()`, `SVDGetConverged()`, `SVDSetWhichSingularTriplets()`
326 @*/
327 35303 PetscErrorCode SVDGetSingularTriplet(SVD svd,PetscInt i,PetscReal *sigma,Vec u,Vec v)
328 {
329 35303 PetscInt M,N;
330 35303 Vec w;
331
332
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
35303 PetscFunctionBegin;
333
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.
35303 PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
334
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.
35303 PetscValidLogicalCollectiveInt(svd,i,2);
335
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
35303 SVDCheckSolved(svd,1);
336
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.
35303 if (u) { PetscValidHeaderSpecific(u,VEC_CLASSID,4); PetscCheckSameComm(svd,1,u,4); }
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.
35303 if (v) { PetscValidHeaderSpecific(v,VEC_CLASSID,5); PetscCheckSameComm(svd,1,v,5); }
338
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
35303 PetscCheck(i>=0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative");
339
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
35303 PetscCheck(i<svd->nconv,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see SVDGetConverged()");
340
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
35303 if (sigma) *sigma = svd->sigma[svd->perm[i]];
341
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
35303 if (u || v) {
342
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
24277 if (!svd->isgeneralized) {
343
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.
18244 PetscCall(MatGetSize(svd->OP,&M,&N));
344
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
18244 if (M<N) { w = u; u = v; v = w; }
345 }
346
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.
24277 PetscCall(SVDComputeVectors(svd));
347
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.
24277 if (u) PetscCall(BVCopyVec(svd->U,svd->perm[i],u));
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.
24277 if (v) PetscCall(BVCopyVec(svd->V,svd->perm[i],v));
349 }
350
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.
7795 PetscFunctionReturn(PETSC_SUCCESS);
351 }
352
353 /*
354 SVDComputeResidualNorms_Standard - Computes the norms of the left and
355 right residuals associated with the i-th computed singular triplet.
356
357 Input Parameters:
358 sigma - singular value
359 u,v - singular vectors
360 x,y - two work vectors with the same dimensions as u,v
361 */
362 5945 static PetscErrorCode SVDComputeResidualNorms_Standard(SVD svd,PetscReal sigma,Vec u,Vec v,Vec x,Vec y,PetscReal *norm1,PetscReal *norm2)
363 {
364 5945 PetscInt M,N;
365
366
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
5945 PetscFunctionBegin;
367 /* norm1 = ||A*v-sigma*u||_2 */
368
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
5945 if (norm1) {
369
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.
5945 PetscCall(MatMult(svd->OP,v,x));
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.
5945 PetscCall(VecAXPY(x,-sigma,u));
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.
5945 PetscCall(VecNorm(x,NORM_2,norm1));
372 }
373 /* norm2 = ||A^T*u-sigma*v||_2 */
374
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
5945 if (norm2) {
375
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.
5945 PetscCall(MatGetSize(svd->OP,&M,&N));
376
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.
5945 if (M<N) PetscCall(MatMult(svd->A,u,y));
377
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.
5218 else PetscCall(MatMult(svd->AT,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.
5945 PetscCall(VecAXPY(y,-sigma,v));
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.
5945 PetscCall(VecNorm(y,NORM_2,norm2));
380 }
381
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.
932 PetscFunctionReturn(PETSC_SUCCESS);
382 }
383
384 /*
385 SVDComputeResidualNorms_Generalized - In GSVD, compute the residual norms
386 norm1 = ||s^2*A'*u-c*B'*B*x||_2 and norm2 = ||c^2*B'*v-s*A'*A*x||_2.
387
388 Input Parameters:
389 sigma - singular value
390 uv - left singular vectors [u;v]
391 x - right singular vector
392 y,z - two work vectors with the same dimension as x
393 */
394 4063 static PetscErrorCode SVDComputeResidualNorms_Generalized(SVD svd,PetscReal sigma,Vec uv,Vec x,Vec y,Vec z,PetscReal *norm1,PetscReal *norm2)
395 {
396 4063 Vec u,v,ax,bx,nest,aux[2];
397 4063 PetscReal c,s;
398
399
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
4063 PetscFunctionBegin;
400
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.
4063 PetscCall(MatCreateVecs(svd->OP,NULL,&u));
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.
4063 PetscCall(MatCreateVecs(svd->OPb,NULL,&v));
402 4063 aux[0] = u;
403 4063 aux[1] = v;
404
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.
4063 PetscCall(VecCreateNest(PetscObjectComm((PetscObject)svd),2,NULL,aux,&nest));
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.
4063 PetscCall(VecCopy(uv,nest));
406
407 4063 s = 1.0/PetscSqrtReal(1.0+sigma*sigma);
408 4063 c = sigma*s;
409
410 /* norm1 = ||s^2*A'*u-c*B'*B*x||_2 */
411
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
4063 if (norm1) {
412
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.
4063 PetscCall(VecDuplicate(v,&bx));
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.
4063 PetscCall(MatMultHermitianTranspose(svd->OP,u,z));
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.
4063 PetscCall(MatMult(svd->OPb,x,bx));
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.
4063 PetscCall(MatMultHermitianTranspose(svd->OPb,bx,y));
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.
4063 PetscCall(VecAXPBY(y,s*s,-c,z));
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.
4063 PetscCall(VecNorm(y,NORM_2,norm1));
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.
4063 PetscCall(VecDestroy(&bx));
419 }
420 /* norm2 = ||c^2*B'*v-s*A'*A*x||_2 */
421
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
4063 if (norm2) {
422
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.
4063 PetscCall(VecDuplicate(u,&ax));
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.
4063 PetscCall(MatMultHermitianTranspose(svd->OPb,v,z));
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.
4063 PetscCall(MatMult(svd->OP,x,ax));
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.
4063 PetscCall(MatMultHermitianTranspose(svd->OP,ax,y));
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.
4063 PetscCall(VecAXPBY(y,c*c,-s,z));
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.
4063 PetscCall(VecNorm(y,NORM_2,norm2));
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.
4063 PetscCall(VecDestroy(&ax));
429 }
430
431
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.
4063 PetscCall(VecDestroy(&v));
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.
4063 PetscCall(VecDestroy(&u));
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.
4063 PetscCall(VecDestroy(&nest));
434
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.
802 PetscFunctionReturn(PETSC_SUCCESS);
435 }
436
437 /*
438 SVDComputeResidualNorms_Hyperbolic - Computes the norms of the left and
439 right residuals associated with the i-th computed singular triplet.
440
441 Input Parameters:
442 sigma - singular value
443 sign - corresponding element of the signature Omega2
444 u,v - singular vectors
445 x,y,z - three work vectors with the same dimensions as u,v,u
446 */
447 4176 static PetscErrorCode SVDComputeResidualNorms_Hyperbolic(SVD svd,PetscReal sigma,PetscReal sign,Vec u,Vec v,Vec x,Vec y,Vec z,PetscReal *norm1,PetscReal *norm2)
448 {
449 4176 PetscInt M,N;
450
451
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
4176 PetscFunctionBegin;
452 /* norm1 = ||A*v-sigma*u||_2 */
453
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
4176 if (norm1) {
454
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.
4176 PetscCall(MatMult(svd->OP,v,x));
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.
4176 PetscCall(VecAXPY(x,-sigma,u));
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.
4176 PetscCall(VecNorm(x,NORM_2,norm1));
457 }
458 /* norm2 = ||A^T*Omega*u-sigma*sign*v||_2 */
459
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
4176 if (norm2) {
460
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.
4176 PetscCall(MatGetSize(svd->OP,&M,&N));
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.
4176 PetscCall(VecPointwiseMult(z,u,svd->omega));
462
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.
4176 if (M<N) PetscCall(MatMult(svd->A,z,y));
463
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.
3888 else PetscCall(MatMult(svd->AT,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.
4176 PetscCall(VecAXPY(y,-sigma*sign,v));
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.
4176 PetscCall(VecNorm(y,NORM_2,norm2));
466 }
467
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.
1276 PetscFunctionReturn(PETSC_SUCCESS);
468 }
469
470 /*@
471 SVDComputeError - Computes the error (based on the residual norm) associated
472 with the `i`-th singular triplet.
473
474 Collective
475
476 Input Parameters:
477 + svd - the singular value solver context
478 . i - the solution index
479 - type - the type of error to compute, see `SVDErrorType`
480
481 Output Parameter:
482 . error - the error
483
484 Notes:
485 The error can be computed in various ways, all of them based on the residual
486 norm obtained as $\sqrt{\eta_1^2+\eta_2^2}$ with $\eta_1 = \|Av-\sigma u\|_2$ and
487 $\eta_2 = \|A^*u-\sigma v\|_2$, where $(\sigma,u,v)$ is the approximate singular
488 triplet.
489
490 In the case of the GSVD, the two components of the residual norm are
491 $\eta_1 = \|s^2 A^*u-cB^*Bx\|_2$ and $\eta_2 = ||c^2 B^*v-sA^*Ax||_2$, where
492 $(\sigma,u,v,x)$ is the approximate generalized singular quadruple, with
493 $\sigma=c/s$.
494
495 Level: beginner
496
497 .seealso: [](ch:svd), `SVDErrorType`, `SVDSolve()`
498 @*/
499 14184 PetscErrorCode SVDComputeError(SVD svd,PetscInt i,SVDErrorType type,PetscReal *error)
500 {
501 14184 PetscReal sigma,norm1,norm2,c,s;
502 14184 Vec u=NULL,v=NULL,x=NULL,y=NULL,z=NULL;
503 14184 PetscReal vecnorm=1.0;
504
505
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
14184 PetscFunctionBegin;
506
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.
14184 PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
507
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.
14184 PetscValidLogicalCollectiveInt(svd,i,2);
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.
14184 PetscValidLogicalCollectiveEnum(svd,type,3);
509
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.
14184 PetscAssertPointer(error,4);
510
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
14184 SVDCheckSolved(svd,1);
511
512 /* allocate work vectors */
513
3/4
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
14184 switch (svd->problem_type) {
514 5945 case SVD_STANDARD:
515
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.
5945 PetscCall(SVDSetWorkVecs(svd,2,2));
516 5945 u = svd->workl[0];
517 5945 v = svd->workr[0];
518 5945 x = svd->workl[1];
519 5945 y = svd->workr[1];
520 5945 break;
521 4063 case SVD_GENERALIZED:
522
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.
4063 PetscCall(SVDSetWorkVecs(svd,1,3));
523 4063 u = svd->workl[0];
524 4063 v = svd->workr[0];
525 4063 x = svd->workr[1];
526 4063 y = svd->workr[2];
527 4063 break;
528 4176 case SVD_HYPERBOLIC:
529
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.
4176 PetscCall(SVDSetWorkVecs(svd,3,2));
530 4176 u = svd->workl[0];
531 4176 v = svd->workr[0];
532 4176 x = svd->workl[1];
533 4176 y = svd->workr[1];
534 4176 z = svd->workl[2];
535 4176 break;
536 }
537
538 /* compute residual norm */
539
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.
14184 PetscCall(SVDGetSingularTriplet(svd,i,&sigma,u,v));
540
3/4
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
14184 switch (svd->problem_type) {
541 5945 case SVD_STANDARD:
542
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.
5945 PetscCall(SVDComputeResidualNorms_Standard(svd,sigma,u,v,x,y,&norm1,&norm2));
543 break;
544 4063 case SVD_GENERALIZED:
545
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.
4063 PetscCall(SVDComputeResidualNorms_Generalized(svd,sigma,u,v,x,y,&norm1,&norm2));
546 break;
547 4176 case SVD_HYPERBOLIC:
548
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.
4176 PetscCall(SVDComputeResidualNorms_Hyperbolic(svd,sigma,svd->sign[svd->perm[i]],u,v,x,y,z,&norm1,&norm2));
549 break;
550 }
551 14184 *error = SlepcAbs(norm1,norm2);
552
553 /* compute 2-norm of eigenvector of the cyclic form */
554
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
14184 if (type!=SVD_ERROR_ABSOLUTE) {
555
3/4
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
14050 switch (svd->problem_type) {
556 5811 case SVD_STANDARD:
557 5811 vecnorm = PETSC_SQRT2;
558 5811 break;
559 4063 case SVD_GENERALIZED:
560
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.
4063 PetscCall(VecNorm(v,NORM_2,&vecnorm));
561 4063 vecnorm = PetscSqrtReal(1.0+vecnorm*vecnorm);
562 4063 break;
563 4176 case SVD_HYPERBOLIC:
564
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.
4176 PetscCall(VecNorm(u,NORM_2,&vecnorm));
565 4176 vecnorm = PetscSqrtReal(1.0+vecnorm*vecnorm);
566 4176 break;
567 }
568 26 }
569
570 /* compute error */
571
3/4
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
14076 switch (type) {
572 case SVD_ERROR_ABSOLUTE:
573 break;
574 5737 case SVD_ERROR_RELATIVE:
575
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
5737 if (svd->isgeneralized) {
576 s = 1.0/PetscSqrtReal(1.0+sigma*sigma);
577 c = sigma*s;
578 norm1 /= c*vecnorm;
579 norm2 /= s*vecnorm;
580 *error = PetscMax(norm1,norm2);
581 5737 } else *error /= sigma*vecnorm;
582 break;
583 8313 case SVD_ERROR_NORM:
584
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.
8313 if (!svd->nrma) PetscCall(MatNorm(svd->OP,NORM_INFINITY,&svd->nrma));
585
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.
8313 if (svd->isgeneralized && !svd->nrmb) PetscCall(MatNorm(svd->OPb,NORM_INFINITY,&svd->nrmb));
586
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
8313 *error /= PetscMax(svd->nrma,svd->nrmb)*vecnorm;
587 8313 break;
588 default:
589 SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
590 }
591
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.
3010 PetscFunctionReturn(PETSC_SUCCESS);
592 }
593