GCC Code Coverage Report


Directory: ./
File: src/svd/impls/lanczos/gklanczos.c
Date: 2026-02-22 03:58:10
Exec Total Coverage
Lines: 249 265 94.0%
Functions: 12 13 92.3%
Branches: 685 1194 57.4%

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 SLEPc singular value solver: "lanczos"
12
13 Method: Explicitly restarted Lanczos
14
15 Algorithm:
16
17 Golub-Kahan-Lanczos bidiagonalization with explicit restart.
18
19 References:
20
21 [1] G.H. Golub and W. Kahan, "Calculating the singular values
22 and pseudo-inverse of a matrix", SIAM J. Numer. Anal. Ser.
23 B 2:205-224, 1965.
24
25 [2] V. Hernandez, J.E. Roman, and A. Tomas, "A robust and
26 efficient parallel SVD solver based on restarted Lanczos
27 bidiagonalization", Elec. Trans. Numer. Anal. 31:68-85,
28 2008.
29 */
30
31 #include <slepc/private/svdimpl.h> /*I "slepcsvd.h" I*/
32
33 typedef struct {
34 PetscBool oneside;
35 } SVD_LANCZOS;
36
37 171 static PetscErrorCode SVDSetUp_Lanczos(SVD svd)
38 {
39 171 SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
40 171 PetscInt N;
41
42
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
171 PetscFunctionBegin;
43
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
171 SVDCheckStandard(svd);
44
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
171 SVDCheckDefinite(svd);
45
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
171 PetscCall(MatGetSize(svd->A,NULL,&N));
46
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
171 PetscCall(SVDSetDimensions_Default(svd));
47
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
171 PetscCheck(svd->ncv<=svd->nsv+svd->mpd,PetscObjectComm((PetscObject)svd),PETSC_ERR_USER_INPUT,"The value of ncv must not be larger than nev+mpd");
48
4/4
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 7 times.
248 if (svd->max_it==PETSC_DETERMINE) svd->max_it = PetscMax(N/svd->ncv,100)*((svd->stop==SVD_STOP_THRESHOLD)?10:1);
49 171 svd->leftbasis = PetscNot(lanczos->oneside);
50
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
171 PetscCall(SVDAllocateSolution(svd,1));
51
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
171 PetscCall(DSSetType(svd->ds,DSSVD));
52
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
171 PetscCall(DSSetCompact(svd->ds,PETSC_TRUE));
53
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
171 PetscCall(DSSetExtraRow(svd->ds,PETSC_TRUE));
54
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
171 PetscCall(DSAllocate(svd->ds,svd->ncv+1));
55
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.
41 PetscFunctionReturn(PETSC_SUCCESS);
56 }
57
58 10165 PetscErrorCode SVDTwoSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt k,PetscInt *n,PetscBool *breakdown)
59 {
60 10165 PetscInt i;
61 10165 Vec u,v;
62 10165 PetscBool lindep=PETSC_FALSE;
63
64
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
10165 PetscFunctionBegin;
65
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
10165 PetscCall(BVGetColumn(svd->V,k,&v));
66
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
10165 PetscCall(BVGetColumn(svd->U,k,&u));
67
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
10165 PetscCall(MatMult(svd->A,v,u));
68
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
10165 PetscCall(BVRestoreColumn(svd->V,k,&v));
69
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
10165 PetscCall(BVRestoreColumn(svd->U,k,&u));
70
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
10165 PetscCall(BVOrthonormalizeColumn(svd->U,k,PETSC_FALSE,alpha+k,&lindep));
71
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
10165 if (PetscUnlikely(lindep)) {
72 *n = k;
73 if (breakdown) *breakdown = lindep;
74 PetscFunctionReturn(PETSC_SUCCESS);
75 }
76
77
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
73345 for (i=k+1;i<*n;i++) {
78
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
63180 PetscCall(BVGetColumn(svd->V,i,&v));
79
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
63180 PetscCall(BVGetColumn(svd->U,i-1,&u));
80
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
63180 PetscCall(MatMult(svd->AT,u,v));
81
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
63180 PetscCall(BVRestoreColumn(svd->V,i,&v));
82
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
63180 PetscCall(BVRestoreColumn(svd->U,i-1,&u));
83
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
63180 PetscCall(BVOrthonormalizeColumn(svd->V,i,PETSC_FALSE,beta+i-1,&lindep));
84
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
63180 if (PetscUnlikely(lindep)) {
85 *n = i;
86 break;
87 }
88
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
63180 PetscCall(BVGetColumn(svd->V,i,&v));
89
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
63180 PetscCall(BVGetColumn(svd->U,i,&u));
90
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
63180 PetscCall(MatMult(svd->A,v,u));
91
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
63180 PetscCall(BVRestoreColumn(svd->V,i,&v));
92
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
63180 PetscCall(BVRestoreColumn(svd->U,i,&u));
93
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
63180 PetscCall(BVOrthonormalizeColumn(svd->U,i,PETSC_FALSE,alpha+i,&lindep));
94
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
63180 if (PetscUnlikely(lindep)) {
95 *n = i;
96 break;
97 }
98 }
99
100
1/2
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
10165 if (!lindep) {
101
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
10165 PetscCall(BVGetColumn(svd->V,*n,&v));
102
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
10165 PetscCall(BVGetColumn(svd->U,*n-1,&u));
103
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
10165 PetscCall(MatMult(svd->AT,u,v));
104
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
10165 PetscCall(BVRestoreColumn(svd->V,*n,&v));
105
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
10165 PetscCall(BVRestoreColumn(svd->U,*n-1,&u));
106
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
10165 PetscCall(BVOrthogonalizeColumn(svd->V,*n,NULL,beta+*n-1,&lindep));
107 }
108
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
10165 if (breakdown) *breakdown = lindep;
109
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.
2477 PetscFunctionReturn(PETSC_SUCCESS);
110 }
111
112 506 static PetscErrorCode SVDOneSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,Vec u,Vec u_1,PetscInt k,PetscInt n,PetscScalar* work)
113 {
114 506 PetscInt i,bvl,bvk;
115 506 PetscReal a,b;
116 506 Vec z,temp;
117
118
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
506 PetscFunctionBegin;
119
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(BVGetActiveColumns(V,&bvl,&bvk));
120
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(BVGetColumn(V,k,&z));
121
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(MatMult(svd->A,z,u));
122
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(BVRestoreColumn(V,k,&z));
123
124
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
4394 for (i=k+1;i<n;i++) {
125
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(BVGetColumn(V,i,&z));
126
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(MatMult(svd->AT,u,z));
127
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(BVRestoreColumn(V,i,&z));
128
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(VecNormBegin(u,NORM_2,&a));
129
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(BVSetActiveColumns(V,0,i));
130
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(BVDotColumnBegin(V,i,work));
131
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(VecNormEnd(u,NORM_2,&a));
132
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(BVDotColumnEnd(V,i,work));
133
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(VecScale(u,1.0/a));
134
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(BVMultColumn(V,-1.0/a,1.0/a,i,work));
135
136 /* h = V^* z, z = z - V h */
137
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(BVDotColumn(V,i,work));
138
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(BVMultColumn(V,-1.0,1.0,i,work));
139
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(BVNormColumn(V,i,NORM_2,&b));
140
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
3888 PetscCheck(PetscAbsReal(b)>10*PETSC_MACHINE_EPSILON,PetscObjectComm((PetscObject)svd),PETSC_ERR_PLIB,"Recurrence generated a zero vector; use a two-sided variant");
141
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(BVScaleColumn(V,i,1.0/b));
142
143
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(BVGetColumn(V,i,&z));
144
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(MatMult(svd->A,z,u_1));
145
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(BVRestoreColumn(V,i,&z));
146
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3888 PetscCall(VecAXPY(u_1,-b,u));
147 3888 alpha[i-1] = a;
148 3888 beta[i-1] = b;
149 3888 temp = u;
150 3888 u = u_1;
151 3888 u_1 = temp;
152 }
153
154
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(BVGetColumn(V,n,&z));
155
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(MatMult(svd->AT,u,z));
156
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(BVRestoreColumn(V,n,&z));
157
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(VecNormBegin(u,NORM_2,&a));
158
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(BVDotColumnBegin(V,n,work));
159
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(VecNormEnd(u,NORM_2,&a));
160
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(BVDotColumnEnd(V,n,work));
161
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(VecScale(u,1.0/a));
162
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(BVMultColumn(V,-1.0/a,1.0/a,n,work));
163
164 /* h = V^* z, z = z - V h */
165
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(BVDotColumn(V,n,work));
166
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(BVMultColumn(V,-1.0,1.0,n,work));
167
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(BVNormColumn(V,i,NORM_2,&b));
168
169 506 alpha[n-1] = a;
170 506 beta[n-1] = b;
171
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(BVSetActiveColumns(V,bvl,bvk));
172
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.
108 PetscFunctionReturn(PETSC_SUCCESS);
173 }
174
175 /*
176 SVDKrylovConvergence - Implements the loop that checks for convergence
177 in Krylov methods.
178
179 Input Parameters:
180 svd - the solver; some error estimates are updated in svd->errest
181 getall - whether all residuals must be computed
182 kini - initial value of k (the loop variable)
183 nits - number of iterations of the loop
184 normr - norm of triangular factor of qr([A;B]), used only in GSVD
185
186 Output Parameter:
187 kout - the first index where the convergence test failed
188 */
189 19901 PetscErrorCode SVDKrylovConvergence(SVD svd,PetscBool getall,PetscInt kini,PetscInt nits,PetscReal normr,PetscInt *kout)
190 {
191 19901 PetscInt k,marker,ld;
192 19901 PetscReal *alpha,*beta,*betah,resnorm;
193 19901 PetscBool extra;
194
195
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
19901 PetscFunctionBegin;
196
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
19901 if (PetscUnlikely(svd->conv == SVD_CONV_MAXIT && svd->its >= svd->max_it)) *kout = svd->nsv;
197 else {
198
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
19901 PetscCall(DSGetLeadingDimension(svd->ds,&ld));
199
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
19901 PetscCall(DSGetExtraRow(svd->ds,&extra));
200
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
19901 PetscCheck(extra,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Only implemented for DS with extra row");
201 19901 marker = -1;
202
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
19901 if (svd->trackall) getall = PETSC_TRUE;
203
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
19901 PetscCall(DSGetArrayReal(svd->ds,DS_MAT_T,&alpha));
204 19901 beta = alpha + ld;
205 19901 betah = alpha + 2*ld;
206
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
25583 for (k=kini;k<kini+nits;k++) {
207
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
25539 if (svd->isgeneralized) resnorm = SlepcAbs(beta[k],betah[k])*normr;
208 16699 else resnorm = PetscAbsReal(beta[k]);
209
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
25539 PetscCall((*svd->converged)(svd,svd->sigma[k],resnorm,&svd->errest[k],svd->convergedctx));
210
3/4
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 8 times.
25539 if (marker==-1 && svd->errest[k] >= svd->tol) marker = k;
211
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
25539 if (marker!=-1 && !getall) break;
212 }
213
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
19901 PetscCall(DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha));
214
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
19901 if (marker!=-1) k = marker;
215 19901 *kout = k;
216 }
217
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.
4957 PetscFunctionReturn(PETSC_SUCCESS);
218 }
219
220 171 static PetscErrorCode SVDSolve_Lanczos(SVD svd)
221 {
222 171 SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
223 171 PetscReal *alpha,*beta;
224 171 PetscScalar *swork,*w,*P,*aux1,*aux2;
225 171 PetscInt i,k,j,nv,ld;
226 171 Vec u=NULL,u_1=NULL;
227 171 Mat U,V;
228
229
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
171 PetscFunctionBegin;
230 /* allocate working space */
231
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
171 PetscCall(DSGetLeadingDimension(svd->ds,&ld));
232
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
171 PetscCall(PetscMalloc2(ld,&w,svd->ncv,&swork));
233
234
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
171 if (lanczos->oneside) {
235
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
18 PetscCall(MatCreateVecs(svd->A,NULL,&u));
236
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
18 PetscCall(MatCreateVecs(svd->A,NULL,&u_1));
237 }
238
239 /* normalize start vector */
240
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
171 if (!svd->nini) {
241
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
135 PetscCall(BVSetRandomColumn(svd->V,0));
242
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
135 PetscCall(BVOrthonormalizeColumn(svd->V,0,PETSC_TRUE,NULL,NULL));
243 }
244
245 5405 while (svd->reason == SVD_CONVERGED_ITERATING) {
246 5234 svd->its++;
247
248 /* inner loop */
249 5234 nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
250
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5234 PetscCall(DSGetArrayReal(svd->ds,DS_MAT_T,&alpha));
251 5234 beta = alpha + ld;
252
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 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.
5234 if (lanczos->oneside) PetscCall(SVDOneSideLanczos(svd,alpha,beta,svd->V,u,u_1,svd->nconv,nv,swork));
253 else {
254
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4728 PetscCall(SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv,&nv,NULL));
255
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4728 PetscCall(BVSetActiveColumns(svd->U,svd->nconv,nv));
256 }
257
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5234 PetscCall(DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha));
258
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5234 PetscCall(BVSetActiveColumns(svd->V,svd->nconv,nv));
259
260 /* compute SVD of bidiagonal matrix */
261
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5234 PetscCall(DSSetDimensions(svd->ds,nv,svd->nconv,0));
262
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5234 PetscCall(DSSVDSetDimensions(svd->ds,nv));
263
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5234 PetscCall(DSSetState(svd->ds,DS_STATE_INTERMEDIATE));
264
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5234 PetscCall(DSSolve(svd->ds,w,NULL));
265
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5234 PetscCall(DSSort(svd->ds,w,NULL,NULL,NULL,NULL));
266
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5234 PetscCall(DSUpdateExtraRow(svd->ds));
267
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5234 PetscCall(DSSynchronize(svd->ds,w,NULL));
268
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
54449 for (i=svd->nconv;i<nv;i++) svd->sigma[i] = PetscRealPart(w[i]);
269
270 /* check convergence */
271
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5234 PetscCall(SVDKrylovConvergence(svd,PETSC_FALSE,svd->nconv,nv-svd->nconv,1.0,&k));
272
4/4
✓ Branch 0 taken 7 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 7 times.
✓ Branch 3 taken 7 times.
5234 SVDSetCtxThreshold(svd,svd->sigma,k);
273
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5234 PetscCall((*svd->stopping)(svd,svd->its,svd->max_it,k,svd->nsv,&svd->reason,svd->stoppingctx));
274
275 /* compute restart vector */
276
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
5234 if (svd->reason == SVD_CONVERGED_ITERATING) {
277
1/2
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
5063 if (k<nv) {
278
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5063 PetscCall(DSGetArray(svd->ds,DS_MAT_V,&P));
279
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
52734 for (j=svd->nconv;j<nv;j++) swork[j-svd->nconv] = PetscConj(P[j+k*ld]);
280
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5063 PetscCall(DSRestoreArray(svd->ds,DS_MAT_V,&P));
281
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5063 PetscCall(BVMultColumn(svd->V,1.0,0.0,nv,swork));
282 } else {
283 /* all approximations have converged, generate a new initial vector */
284 PetscCall(BVSetRandomColumn(svd->V,nv));
285 PetscCall(BVOrthonormalizeColumn(svd->V,nv,PETSC_FALSE,NULL,NULL));
286 }
287 }
288
289 /* compute converged singular vectors */
290
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5234 PetscCall(DSGetMat(svd->ds,DS_MAT_V,&V));
291
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5234 PetscCall(BVMultInPlace(svd->V,V,svd->nconv,k));
292
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5234 PetscCall(DSRestoreMat(svd->ds,DS_MAT_V,&V));
293
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
5234 if (!lanczos->oneside) {
294
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4728 PetscCall(DSGetMat(svd->ds,DS_MAT_U,&U));
295
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4728 PetscCall(BVMultInPlace(svd->U,U,svd->nconv,k));
296
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4728 PetscCall(DSRestoreMat(svd->ds,DS_MAT_U,&U));
297 }
298
299
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
5234 if (svd->reason == SVD_CONVERGED_ITERATING) {
300
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5063 PetscCall(BVCopyColumn(svd->V,nv,k)); /* copy restart vector from the last column */
301
4/4
✓ Branch 0 taken 7 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 7 times.
✓ Branch 3 taken 7 times.
5063 if (svd->stop==SVD_STOP_THRESHOLD && nv-k<5) { /* reallocate */
302 7 svd->ncv = svd->mpd+k;
303
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7 PetscCall(SVDReallocateSolution(svd,svd->ncv+1));
304
2/2
✓ Branch 0 taken 7 times.
✓ Branch 1 taken 7 times.
57 for (i=nv;i<svd->ncv;i++) svd->perm[i] = i;
305
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7 PetscCall(DSReallocate(svd->ds,svd->ncv+1));
306 7 aux1 = w;
307 7 aux2 = swork;
308
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7 PetscCall(PetscMalloc2(svd->ncv+1,&w,svd->ncv,&swork));
309
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7 PetscCall(PetscArraycpy(w,aux1,ld));
310
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7 PetscCall(PetscFree2(aux1,aux2));
311
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7 PetscCall(DSGetLeadingDimension(svd->ds,&ld));
312 }
313 }
314
315 5234 svd->nconv = k;
316
7/8
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 6 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
5405 PetscCall(SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv));
317 }
318
319 /* free working space */
320
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
171 PetscCall(VecDestroy(&u));
321
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
171 PetscCall(VecDestroy(&u_1));
322
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
171 PetscCall(PetscFree2(w,swork));
323
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.
41 PetscFunctionReturn(PETSC_SUCCESS);
324 }
325
326 123 static PetscErrorCode SVDSetFromOptions_Lanczos(SVD svd,PetscOptionItems PetscOptionsObject)
327 {
328 123 PetscBool set,val;
329 123 SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
330
331
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
123 PetscFunctionBegin;
332
1/12
✗ Branch 0 not taken.
✓ Branch 1 taken 8 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.
123 PetscOptionsHeadBegin(PetscOptionsObject,"SVD Lanczos Options");
333
334
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
123 PetscCall(PetscOptionsBool("-svd_lanczos_oneside","Use one-side reorthogonalization","SVDLanczosSetOneSide",lanczos->oneside,&val,&set));
335
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 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.
123 if (set) PetscCall(SVDLanczosSetOneSide(svd,val));
336
337
1/14
✗ Branch 0 not taken.
✓ 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.
123 PetscOptionsHeadEnd();
338
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.
29 PetscFunctionReturn(PETSC_SUCCESS);
339 }
340
341 18 static PetscErrorCode SVDLanczosSetOneSide_Lanczos(SVD svd,PetscBool oneside)
342 {
343 18 SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
344
345
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
18 PetscFunctionBegin;
346
1/2
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
18 if (lanczos->oneside != oneside) {
347 18 lanczos->oneside = oneside;
348 18 svd->state = SVD_STATE_INITIAL;
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.
18 PetscFunctionReturn(PETSC_SUCCESS);
351 }
352
353 /*@
354 SVDLanczosSetOneSide - Indicate if the variant of the Lanczos method
355 to be used is one-sided or two-sided.
356
357 Logically Collective
358
359 Input Parameters:
360 + svd - the singular value solver context
361 - oneside - boolean flag indicating if the method is one-sided or not
362
363 Options Database Key:
364 . -svd_lanczos_oneside - enable the one-sided variant
365
366 Note:
367 By default, a two-sided variant is selected, which is sometimes slightly
368 more robust. However, the one-sided variant is faster because it avoids
369 the orthogonalization associated to left singular vectors. It also saves
370 the memory required for storing such vectors. See more details in
371 {cite:p}`Her07c,Her08`.
372
373 Level: advanced
374
375 .seealso: [](ch:svd), `SVDTRLANCZOS`, `SVDTRLanczosSetOneSide()`
376 @*/
377 18 PetscErrorCode SVDLanczosSetOneSide(SVD svd,PetscBool oneside)
378 {
379
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
18 PetscFunctionBegin;
380
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.
18 PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
381
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.
18 PetscValidLogicalCollectiveBool(svd,oneside,2);
382
8/14
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 8 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 8 times.
✓ 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.
18 PetscTryMethod(svd,"SVDLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
383
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.
18 PetscFunctionReturn(PETSC_SUCCESS);
384 }
385
386 16 static PetscErrorCode SVDLanczosGetOneSide_Lanczos(SVD svd,PetscBool *oneside)
387 {
388 16 SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
389
390
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
16 PetscFunctionBegin;
391 16 *oneside = lanczos->oneside;
392
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.
16 PetscFunctionReturn(PETSC_SUCCESS);
393 }
394
395 /*@
396 SVDLanczosGetOneSide - Gets if the variant of the Lanczos method
397 to be used is one-sided or two-sided.
398
399 Not Collective
400
401 Input Parameter:
402 . svd - the singular value solver context
403
404 Output Parameter:
405 . oneside - boolean flag indicating if the method is one-sided or not
406
407 Level: advanced
408
409 .seealso: [](ch:svd), `SVDTRLANCZOS`, `SVDLanczosSetOneSide()`
410 @*/
411 16 PetscErrorCode SVDLanczosGetOneSide(SVD svd,PetscBool *oneside)
412 {
413
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
16 PetscFunctionBegin;
414
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.
16 PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
415
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.
16 PetscAssertPointer(oneside,2);
416
9/16
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 8 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
16 PetscUseMethod(svd,"SVDLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
417
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.
16 PetscFunctionReturn(PETSC_SUCCESS);
418 }
419
420 139 static PetscErrorCode SVDDestroy_Lanczos(SVD svd)
421 {
422
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
139 PetscFunctionBegin;
423
5/8
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
✓ 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.
139 PetscCall(PetscFree(svd->data));
424
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
139 PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",NULL));
425
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
139 PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",NULL));
426
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.
33 PetscFunctionReturn(PETSC_SUCCESS);
427 }
428
429 static PetscErrorCode SVDView_Lanczos(SVD svd,PetscViewer viewer)
430 {
431 SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
432 PetscBool isascii;
433
434 PetscFunctionBegin;
435 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
436 if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer," %s-sided reorthogonalization\n",lanczos->oneside? "one": "two"));
437 PetscFunctionReturn(PETSC_SUCCESS);
438 }
439
440 /*MC
441 SVDLANCZOS - SVDLANCZOS = "lanczos" - A basic Golub-Kahan-Lanczos
442 bidiagonalization method with explicit restart.
443
444 Notes:
445 Only available for standard SVD problems.
446
447 This solver is very basic and is not recommended in general, since it
448 will not be competitive with respect to other solvers.
449
450 The implemented method is Lanczos bidiagonalization with explicit restart
451 and deflation. Generally, it is much better to use Lanczos with implicit
452 restart (also known as thick-restart Lanczos) as implemented in `SVDTRLANCZOS`.
453
454 The implementation includes a one-sided orthogonalization option, and
455 efficient parallel orthogonalization, see the details in {cite:p}`Her07c,Her08`.
456
457 Level: beginner
458
459 .seealso: [](ch:svd), [](#sec:svdback), `SVD`, `SVDType`, `SVDSetType()`, `SVDSetProblemType()`, `SVDTRLANCZOS`
460 M*/
461 139 SLEPC_EXTERN PetscErrorCode SVDCreate_Lanczos(SVD svd)
462 {
463 139 SVD_LANCZOS *ctx;
464
465
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
139 PetscFunctionBegin;
466
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
139 PetscCall(PetscNew(&ctx));
467 139 svd->data = (void*)ctx;
468
469 139 svd->ops->setup = SVDSetUp_Lanczos;
470 139 svd->ops->solve = SVDSolve_Lanczos;
471 139 svd->ops->destroy = SVDDestroy_Lanczos;
472 139 svd->ops->setfromoptions = SVDSetFromOptions_Lanczos;
473 139 svd->ops->view = SVDView_Lanczos;
474 139 svd->ops->computevectors = SVDComputeVectors_Left;
475
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
139 PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",SVDLanczosSetOneSide_Lanczos));
476
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
139 PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",SVDLanczosGetOneSide_Lanczos));
477
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.
33 PetscFunctionReturn(PETSC_SUCCESS);
478 }
479