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 - singular value solver context obtained from SVDCreate() | ||
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, 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 | All the command-line options listed above admit an optional argument specifying | ||
107 | the viewer type and options. For instance, use '-svd_view_mat0 binary:amatrix.bin' | ||
108 | to save the A matrix to a binary file, '-svd_view_values draw' to draw the computed | ||
109 | singular values graphically, or '-svd_error_relative :myerr.m:ascii_matlab' to save | ||
110 | the errors in a file that can be executed in Matlab. | ||
111 | |||
112 | Level: beginner | ||
113 | |||
114 | .seealso: SVDCreate(), SVDSetUp(), SVDDestroy() | ||
115 | @*/ | ||
116 | 2773 | PetscErrorCode SVDSolve(SVD svd) | |
117 | { | ||
118 | 2773 | PetscInt i,m,n,*workperm; | |
119 | |||
120 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
2773 | PetscFunctionBegin; |
121 |
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); |
122 |
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); |
123 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2773 | PetscCall(PetscLogEventBegin(SVD_Solve,svd,0,0,0)); |
124 | |||
125 | /* call setup */ | ||
126 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2773 | PetscCall(SVDSetUp(svd)); |
127 | |||
128 | /* safeguard for matrices with zero rows or columns */ | ||
129 |
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)); |
130 |
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) { |
131 | ✗ | svd->nconv = 0; | |
132 | ✗ | svd->reason = SVD_CONVERGED_TOL; | |
133 | ✗ | svd->state = SVD_STATE_SOLVED; | |
134 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
135 | } | ||
136 | |||
137 | 2773 | svd->its = 0; | |
138 | 2773 | svd->nconv = 0; | |
139 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
61325 | for (i=0;i<svd->ncv;i++) { |
140 | 58552 | svd->sigma[i] = 0.0; | |
141 | 58552 | svd->errest[i] = 0.0; | |
142 | 58552 | svd->perm[i] = i; | |
143 | } | ||
144 |
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")); |
145 | |||
146 |
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) { |
147 | 1526 | case SVD_STANDARD: | |
148 |
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); |
149 | break; | ||
150 | 985 | case SVD_GENERALIZED: | |
151 |
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); |
152 | break; | ||
153 | 262 | case SVD_HYPERBOLIC: | |
154 |
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); |
155 | break; | ||
156 | } | ||
157 | 2773 | svd->state = SVD_STATE_SOLVED; | |
158 | |||
159 | /* sort singular triplets */ | ||
160 |
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)); |
161 | else { | ||
162 |
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)); |
163 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
35568 | for (i=0;i<svd->nconv;i++) workperm[i] = i; |
164 |
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)); |
165 |
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]; |
166 |
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)); |
167 | } | ||
168 |
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)); |
169 | |||
170 | /* various viewers */ | ||
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.
|
2773 | PetscCall(SVDViewFromOptions(svd,NULL,"-svd_view")); |
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.
|
2773 | PetscCall(SVDConvergedReasonViewFromOptions(svd)); |
173 |
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)); |
174 |
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)); |
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(SVDVectorsViewFromOptions(svd)); |
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.
|
2773 | PetscCall(MatViewFromOptions(svd->OP,(PetscObject)svd,"-svd_view_mat0")); |
177 |
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")); |
178 |
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")); |
179 | |||
180 | /* Remove the initial subspaces */ | ||
181 | 2773 | svd->nini = 0; | |
182 | 2773 | svd->ninil = 0; | |
183 |
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); |
184 | } | ||
185 | |||
186 | /*@ | ||
187 | SVDGetIterationNumber - Gets the current iteration number. If the | ||
188 | call to SVDSolve() is complete, then it returns the number of iterations | ||
189 | carried out by the solution method. | ||
190 | |||
191 | Not Collective | ||
192 | |||
193 | Input Parameter: | ||
194 | . svd - the singular value solver context | ||
195 | |||
196 | Output Parameter: | ||
197 | . its - number of iterations | ||
198 | |||
199 | Note: | ||
200 | During the i-th iteration this call returns i-1. If SVDSolve() is | ||
201 | complete, then parameter "its" contains either the iteration number at | ||
202 | which convergence was successfully reached, or failure was detected. | ||
203 | Call SVDGetConvergedReason() to determine if the solver converged or | ||
204 | failed and why. | ||
205 | |||
206 | Level: intermediate | ||
207 | |||
208 | .seealso: SVDGetConvergedReason(), SVDSetTolerances() | ||
209 | @*/ | ||
210 | 729 | PetscErrorCode SVDGetIterationNumber(SVD svd,PetscInt *its) | |
211 | { | ||
212 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
729 | PetscFunctionBegin; |
213 |
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); |
214 |
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); |
215 | 729 | *its = svd->its; | |
216 |
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); |
217 | } | ||
218 | |||
219 | /*@ | ||
220 | SVDGetConvergedReason - Gets the reason why the SVDSolve() iteration was | ||
221 | stopped. | ||
222 | |||
223 | Not Collective | ||
224 | |||
225 | Input Parameter: | ||
226 | . svd - the singular value solver context | ||
227 | |||
228 | Output Parameter: | ||
229 | . reason - negative value indicates diverged, positive value converged | ||
230 | (see SVDConvergedReason) | ||
231 | |||
232 | Options Database Key: | ||
233 | . -svd_converged_reason - print the reason to a viewer | ||
234 | |||
235 | Notes: | ||
236 | Possible values for reason are | ||
237 | + SVD_CONVERGED_TOL - converged up to tolerance | ||
238 | . SVD_CONVERGED_USER - converged due to a user-defined condition | ||
239 | . SVD_CONVERGED_MAXIT - reached the maximum number of iterations with SVD_CONV_MAXIT criterion | ||
240 | . SVD_DIVERGED_ITS - required more than max_it iterations to reach convergence | ||
241 | . SVD_DIVERGED_BREAKDOWN - generic breakdown in method | ||
242 | - SVD_DIVERGED_SYMMETRY_LOST - underlying indefinite eigensolver was not able to keep symmetry | ||
243 | |||
244 | Can only be called after the call to SVDSolve() is complete. | ||
245 | |||
246 | Level: intermediate | ||
247 | |||
248 | .seealso: 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 | Note: | ||
272 | This function should be called after SVDSolve() has finished. | ||
273 | |||
274 | Level: beginner | ||
275 | |||
276 | .seealso: SVDSetDimensions(), SVDSolve(), SVDGetSingularTriplet() | ||
277 | @*/ | ||
278 | 1213 | PetscErrorCode SVDGetConverged(SVD svd,PetscInt *nconv) | |
279 | { | ||
280 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1213 | PetscFunctionBegin; |
281 |
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); |
282 |
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); |
283 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
1213 | SVDCheckSolved(svd,1); |
284 | 1213 | *nconv = svd->nconv; | |
285 |
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); |
286 | } | ||
287 | |||
288 | /*@ | ||
289 | SVDGetSingularTriplet - Gets the i-th triplet of the singular value decomposition | ||
290 | as computed by SVDSolve(). The solution consists in the singular value and its left | ||
291 | and right singular vectors. | ||
292 | |||
293 | Collective | ||
294 | |||
295 | Input Parameters: | ||
296 | + svd - singular value solver context | ||
297 | - i - index of the solution | ||
298 | |||
299 | Output Parameters: | ||
300 | + sigma - singular value | ||
301 | . u - left singular vector | ||
302 | - v - right singular vector | ||
303 | |||
304 | Note: | ||
305 | Both u or v can be NULL if singular vectors are not required. | ||
306 | Otherwise, the caller must provide valid Vec objects, i.e., | ||
307 | they must be created by the calling program with e.g. MatCreateVecs(). | ||
308 | |||
309 | The index i should be a value between 0 and nconv-1 (see SVDGetConverged()). | ||
310 | Singular triplets are indexed according to the ordering criterion established | ||
311 | with SVDSetWhichSingularTriplets(). | ||
312 | |||
313 | In the case of GSVD, the solution consists in three vectors u,v,x that are | ||
314 | returned as follows. Vector x is returned in the right singular vector | ||
315 | (argument v) and has length equal to the number of columns of A and B. | ||
316 | The other two vectors are returned stacked on top of each other [u;v] in | ||
317 | the left singular vector argument, with length equal to m+n (number of rows | ||
318 | of A plus number of rows of B). | ||
319 | |||
320 | Level: beginner | ||
321 | |||
322 | .seealso: SVDSolve(), SVDGetConverged(), SVDSetWhichSingularTriplets() | ||
323 | @*/ | ||
324 | 35303 | PetscErrorCode SVDGetSingularTriplet(SVD svd,PetscInt i,PetscReal *sigma,Vec u,Vec v) | |
325 | { | ||
326 | 35303 | PetscInt M,N; | |
327 | 35303 | Vec w; | |
328 | |||
329 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
35303 | PetscFunctionBegin; |
330 |
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); |
331 |
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); |
332 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
35303 | SVDCheckSolved(svd,1); |
333 |
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); } |
334 |
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); } |
335 |
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"); |
336 |
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()"); |
337 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
35303 | if (sigma) *sigma = svd->sigma[svd->perm[i]]; |
338 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
35303 | if (u || v) { |
339 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
24277 | if (!svd->isgeneralized) { |
340 |
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)); |
341 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
18244 | if (M<N) { w = u; u = v; v = w; } |
342 | } | ||
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.
|
24277 | PetscCall(SVDComputeVectors(svd)); |
344 |
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)); |
345 |
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)); |
346 | } | ||
347 |
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); |
348 | } | ||
349 | |||
350 | /* | ||
351 | SVDComputeResidualNorms_Standard - Computes the norms of the left and | ||
352 | right residuals associated with the i-th computed singular triplet. | ||
353 | |||
354 | Input Parameters: | ||
355 | sigma - singular value | ||
356 | u,v - singular vectors | ||
357 | x,y - two work vectors with the same dimensions as u,v | ||
358 | */ | ||
359 | 5945 | static PetscErrorCode SVDComputeResidualNorms_Standard(SVD svd,PetscReal sigma,Vec u,Vec v,Vec x,Vec y,PetscReal *norm1,PetscReal *norm2) | |
360 | { | ||
361 | 5945 | PetscInt M,N; | |
362 | |||
363 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
5945 | PetscFunctionBegin; |
364 | /* norm1 = ||A*v-sigma*u||_2 */ | ||
365 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
5945 | if (norm1) { |
366 |
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)); |
367 |
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)); |
368 |
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)); |
369 | } | ||
370 | /* norm2 = ||A^T*u-sigma*v||_2 */ | ||
371 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
5945 | if (norm2) { |
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.
|
5945 | PetscCall(MatGetSize(svd->OP,&M,&N)); |
373 |
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)); |
374 |
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)); |
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(VecAXPY(y,-sigma,v)); |
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.
|
5945 | PetscCall(VecNorm(y,NORM_2,norm2)); |
377 | } | ||
378 |
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); |
379 | } | ||
380 | |||
381 | /* | ||
382 | SVDComputeResidualNorms_Generalized - In GSVD, compute the residual norms | ||
383 | norm1 = ||s^2*A'*u-c*B'*B*x||_2 and norm2 = ||c^2*B'*v-s*A'*A*x||_2. | ||
384 | |||
385 | Input Parameters: | ||
386 | sigma - singular value | ||
387 | uv - left singular vectors [u;v] | ||
388 | x - right singular vector | ||
389 | y,z - two work vectors with the same dimension as x | ||
390 | */ | ||
391 | 4063 | static PetscErrorCode SVDComputeResidualNorms_Generalized(SVD svd,PetscReal sigma,Vec uv,Vec x,Vec y,Vec z,PetscReal *norm1,PetscReal *norm2) | |
392 | { | ||
393 | 4063 | Vec u,v,ax,bx,nest,aux[2]; | |
394 | 4063 | PetscReal c,s; | |
395 | |||
396 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
4063 | PetscFunctionBegin; |
397 |
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)); |
398 |
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)); |
399 | 4063 | aux[0] = u; | |
400 | 4063 | aux[1] = v; | |
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(VecCreateNest(PetscObjectComm((PetscObject)svd),2,NULL,aux,&nest)); |
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.
|
4063 | PetscCall(VecCopy(uv,nest)); |
403 | |||
404 | 4063 | s = 1.0/PetscSqrtReal(1.0+sigma*sigma); | |
405 | 4063 | c = sigma*s; | |
406 | |||
407 | /* norm1 = ||s^2*A'*u-c*B'*B*x||_2 */ | ||
408 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
4063 | if (norm1) { |
409 |
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)); |
410 |
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)); |
411 |
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)); |
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(MatMultHermitianTranspose(svd->OPb,bx,y)); |
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(VecAXPBY(y,s*s,-c,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(VecNorm(y,NORM_2,norm1)); |
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(VecDestroy(&bx)); |
416 | } | ||
417 | /* norm2 = ||c^2*B'*v-s*A'*A*x||_2 */ | ||
418 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
4063 | if (norm2) { |
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.
|
4063 | PetscCall(VecDuplicate(u,&ax)); |
420 |
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)); |
421 |
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)); |
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(MatMultHermitianTranspose(svd->OP,ax,y)); |
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(VecAXPBY(y,c*c,-s,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(VecNorm(y,NORM_2,norm2)); |
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(VecDestroy(&ax)); |
426 | } | ||
427 | |||
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(&v)); |
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.
|
4063 | PetscCall(VecDestroy(&u)); |
430 |
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)); |
431 |
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); |
432 | } | ||
433 | |||
434 | /* | ||
435 | SVDComputeResidualNorms_Hyperbolic - Computes the norms of the left and | ||
436 | right residuals associated with the i-th computed singular triplet. | ||
437 | |||
438 | Input Parameters: | ||
439 | sigma - singular value | ||
440 | sign - corresponding element of the signature Omega2 | ||
441 | u,v - singular vectors | ||
442 | x,y,z - three work vectors with the same dimensions as u,v,u | ||
443 | */ | ||
444 | 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) | |
445 | { | ||
446 | 4176 | PetscInt M,N; | |
447 | |||
448 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
4176 | PetscFunctionBegin; |
449 | /* norm1 = ||A*v-sigma*u||_2 */ | ||
450 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
4176 | if (norm1) { |
451 |
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)); |
452 |
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)); |
453 |
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)); |
454 | } | ||
455 | /* norm2 = ||A^T*Omega*u-sigma*sign*v||_2 */ | ||
456 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
4176 | if (norm2) { |
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.
|
4176 | PetscCall(MatGetSize(svd->OP,&M,&N)); |
458 |
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)); |
459 |
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)); |
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.
|
3888 | else PetscCall(MatMult(svd->AT,z,y)); |
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(VecAXPY(y,-sigma*sign,v)); |
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.
|
4176 | PetscCall(VecNorm(y,NORM_2,norm2)); |
463 | } | ||
464 |
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); |
465 | } | ||
466 | |||
467 | /*@ | ||
468 | SVDComputeError - Computes the error (based on the residual norm) associated | ||
469 | with the i-th singular triplet. | ||
470 | |||
471 | Collective | ||
472 | |||
473 | Input Parameters: | ||
474 | + svd - the singular value solver context | ||
475 | . i - the solution index | ||
476 | - type - the type of error to compute | ||
477 | |||
478 | Output Parameter: | ||
479 | . error - the error | ||
480 | |||
481 | Notes: | ||
482 | The error can be computed in various ways, all of them based on the residual | ||
483 | norm obtained as sqrt(n1^2+n2^2) with n1 = ||A*v-sigma*u||_2 and | ||
484 | n2 = ||A^T*u-sigma*v||_2, where sigma is the singular value, u is the left | ||
485 | singular vector and v is the right singular vector. | ||
486 | |||
487 | In the case of the GSVD, the two components of the residual norm are | ||
488 | n1 = ||s^2*A'*u-c*B'*B*x||_2 and n2 = ||c^2*B'*v-s*A'*A*x||_2, where [u;v] | ||
489 | are the left singular vectors and x is the right singular vector, with | ||
490 | sigma=c/s. | ||
491 | |||
492 | Level: beginner | ||
493 | |||
494 | .seealso: SVDErrorType, SVDSolve() | ||
495 | @*/ | ||
496 | 14184 | PetscErrorCode SVDComputeError(SVD svd,PetscInt i,SVDErrorType type,PetscReal *error) | |
497 | { | ||
498 | 14184 | PetscReal sigma,norm1,norm2,c,s; | |
499 | 14184 | Vec u=NULL,v=NULL,x=NULL,y=NULL,z=NULL; | |
500 | 14184 | PetscReal vecnorm=1.0; | |
501 | |||
502 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
14184 | PetscFunctionBegin; |
503 |
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); |
504 |
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); |
505 |
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); |
506 |
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); |
507 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
14184 | SVDCheckSolved(svd,1); |
508 | |||
509 | /* allocate work vectors */ | ||
510 |
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) { |
511 | 5945 | case SVD_STANDARD: | |
512 |
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)); |
513 | 5945 | u = svd->workl[0]; | |
514 | 5945 | v = svd->workr[0]; | |
515 | 5945 | x = svd->workl[1]; | |
516 | 5945 | y = svd->workr[1]; | |
517 | 5945 | break; | |
518 | 4063 | case SVD_GENERALIZED: | |
519 |
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)); |
520 | 4063 | u = svd->workl[0]; | |
521 | 4063 | v = svd->workr[0]; | |
522 | 4063 | x = svd->workr[1]; | |
523 | 4063 | y = svd->workr[2]; | |
524 | 4063 | break; | |
525 | 4176 | case SVD_HYPERBOLIC: | |
526 |
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)); |
527 | 4176 | u = svd->workl[0]; | |
528 | 4176 | v = svd->workr[0]; | |
529 | 4176 | x = svd->workl[1]; | |
530 | 4176 | y = svd->workr[1]; | |
531 | 4176 | z = svd->workl[2]; | |
532 | 4176 | break; | |
533 | } | ||
534 | |||
535 | /* compute residual norm */ | ||
536 |
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)); |
537 |
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) { |
538 | 5945 | case SVD_STANDARD: | |
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.
|
5945 | PetscCall(SVDComputeResidualNorms_Standard(svd,sigma,u,v,x,y,&norm1,&norm2)); |
540 | break; | ||
541 | 4063 | case SVD_GENERALIZED: | |
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.
|
4063 | PetscCall(SVDComputeResidualNorms_Generalized(svd,sigma,u,v,x,y,&norm1,&norm2)); |
543 | break; | ||
544 | 4176 | case SVD_HYPERBOLIC: | |
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.
|
4176 | PetscCall(SVDComputeResidualNorms_Hyperbolic(svd,sigma,svd->sign[svd->perm[i]],u,v,x,y,z,&norm1,&norm2)); |
546 | break; | ||
547 | } | ||
548 | 14184 | *error = SlepcAbs(norm1,norm2); | |
549 | |||
550 | /* compute 2-norm of eigenvector of the cyclic form */ | ||
551 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
14184 | if (type!=SVD_ERROR_ABSOLUTE) { |
552 |
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) { |
553 | 5811 | case SVD_STANDARD: | |
554 | 5811 | vecnorm = PETSC_SQRT2; | |
555 | 5811 | break; | |
556 | 4063 | case SVD_GENERALIZED: | |
557 |
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)); |
558 | 4063 | vecnorm = PetscSqrtReal(1.0+vecnorm*vecnorm); | |
559 | 4063 | break; | |
560 | 4176 | case SVD_HYPERBOLIC: | |
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.
|
4176 | PetscCall(VecNorm(u,NORM_2,&vecnorm)); |
562 | 4176 | vecnorm = PetscSqrtReal(1.0+vecnorm*vecnorm); | |
563 | 4176 | break; | |
564 | } | ||
565 | 26 | } | |
566 | |||
567 | /* compute error */ | ||
568 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
|
14076 | switch (type) { |
569 | case SVD_ERROR_ABSOLUTE: | ||
570 | break; | ||
571 | 5737 | case SVD_ERROR_RELATIVE: | |
572 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
5737 | if (svd->isgeneralized) { |
573 | ✗ | s = 1.0/PetscSqrtReal(1.0+sigma*sigma); | |
574 | ✗ | c = sigma*s; | |
575 | ✗ | norm1 /= c*vecnorm; | |
576 | ✗ | norm2 /= s*vecnorm; | |
577 | ✗ | *error = PetscMax(norm1,norm2); | |
578 | 5737 | } else *error /= sigma*vecnorm; | |
579 | break; | ||
580 | 8313 | case SVD_ERROR_NORM: | |
581 |
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)); |
582 |
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)); |
583 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
8313 | *error /= PetscMax(svd->nrma,svd->nrmb)*vecnorm; |
584 | 8313 | break; | |
585 | ✗ | default: | |
586 | ✗ | SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type"); | |
587 | } | ||
588 |
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); |
589 | } | ||
590 |