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 | #include <slepc/private/dsimpl.h> /*I "slepcds.h" I*/ | ||
11 | #include <slepcblaslapack.h> | ||
12 | |||
13 | typedef struct { | ||
14 | PetscInt m; /* number of columns */ | ||
15 | PetscInt p; /* number of rows of B */ | ||
16 | PetscInt tm; /* number of rows of X after truncating */ | ||
17 | PetscInt tp; /* number of rows of V after truncating */ | ||
18 | } DS_GSVD; | ||
19 | |||
20 | 684 | static PetscErrorCode DSAllocate_GSVD(DS ds,PetscInt ld) | |
21 | { | ||
22 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
684 | PetscFunctionBegin; |
23 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
684 | PetscCall(DSAllocateMat_Private(ds,DS_MAT_A)); |
24 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
684 | PetscCall(DSAllocateMat_Private(ds,DS_MAT_B)); |
25 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
684 | PetscCall(DSAllocateMat_Private(ds,DS_MAT_X)); |
26 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
684 | PetscCall(DSAllocateMat_Private(ds,DS_MAT_U)); |
27 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
684 | PetscCall(DSAllocateMat_Private(ds,DS_MAT_V)); |
28 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
684 | PetscCall(DSAllocateMat_Private(ds,DS_MAT_T)); |
29 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
684 | PetscCall(DSAllocateMat_Private(ds,DS_MAT_D)); |
30 |
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.
|
684 | PetscCall(PetscFree(ds->perm)); |
31 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
684 | PetscCall(PetscMalloc1(ld,&ds->perm)); |
32 |
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.
|
140 | PetscFunctionReturn(PETSC_SUCCESS); |
33 | } | ||
34 | |||
35 | /* | ||
36 | In compact form, A is either in form (a) or (b): | ||
37 | |||
38 | (a) (b) | ||
39 | lower bidiagonal with upper arrow (n=m+1) square upper bidiagonal with upper arrow (n=m) | ||
40 | 0 l k m-1 | ||
41 | ----------------------------------------- 0 l k m-1 | ||
42 | |* . | ----------------------------------------- | ||
43 | | * . | |* . | | ||
44 | | * . | | * . | | ||
45 | | * . | | * . | | ||
46 | l |. . . . o o | l |. . . o o | | ||
47 | | o o | | o o | | ||
48 | | o o | | o o | | ||
49 | | o o | | o o | | ||
50 | | o o | | o o | | ||
51 | | o o | | o o | | ||
52 | k |. . . . . . . . . . o | k |. . . . . . . . . o x | | ||
53 | | x x | | x x | | ||
54 | | x x | | x x | | ||
55 | | x x | | x x | | ||
56 | | x x | | x x | | ||
57 | | x x | | x x | | ||
58 | | x x | | x x | | ||
59 | | x x | | x x | | ||
60 | | x x | | x x | | ||
61 | | x x| | x x| | ||
62 | n-1 | x| n-1 | x| | ||
63 | ----------------------------------------- ----------------------------------------- | ||
64 | |||
65 | and B is square bidiagonal with upper arrow (p=m) | ||
66 | |||
67 | 0 l k m-1 | ||
68 | ----------------------------------------- | ||
69 | |* . | | ||
70 | | * . | | ||
71 | | * . | | ||
72 | | * . | | ||
73 | l |. . . . o o | | ||
74 | | o o | | ||
75 | | o o | | ||
76 | | o o | | ||
77 | | o o | | ||
78 | | o o | | ||
79 | k |. . . . . . . . . . o x | | ||
80 | | x x | | ||
81 | | x x | | ||
82 | | x x | | ||
83 | | x x | | ||
84 | | x x | | ||
85 | | x x | | ||
86 | | x x | | ||
87 | | x x| | ||
88 | p-1 | x| | ||
89 | ---------------------------------------- | ||
90 | */ | ||
91 | 300 | static PetscErrorCode DSView_GSVD(DS ds,PetscViewer viewer) | |
92 | { | ||
93 | 300 | DS_GSVD *ctx = (DS_GSVD*)ds->data; | |
94 | 300 | PetscViewerFormat format; | |
95 | 300 | PetscInt i,j,r,k=ds->k,n=ds->n,m=ctx->m,p=ctx->p,rowsa,rowsb,colsa,colsb; | |
96 | 300 | PetscReal *T,*S,value; | |
97 | |||
98 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
300 | PetscFunctionBegin; |
99 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
300 | PetscCall(PetscViewerGetFormat(viewer,&format)); |
100 |
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.
|
300 | if (format == PETSC_VIEWER_ASCII_INFO) PetscFunctionReturn(PETSC_SUCCESS); |
101 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
300 | if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { |
102 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
170 | PetscCall(PetscViewerASCIIPrintf(viewer,"number of columns: %" PetscInt_FMT "\n",m)); |
103 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
170 | PetscCall(PetscViewerASCIIPrintf(viewer,"number of rows of B: %" PetscInt_FMT "\n",p)); |
104 |
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.
|
34 | PetscFunctionReturn(PETSC_SUCCESS); |
105 | } | ||
106 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
130 | PetscCheck(ctx->m,PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the other dimensions with DSGSVDSetDimensions()"); |
107 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
130 | if (ds->compact) { |
108 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T)); |
109 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(DSGetArrayReal(ds,DS_MAT_D,&S)); |
110 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); |
111 | 10 | rowsa = n; | |
112 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
10 | colsa = ds->extrarow? m+1: m; |
113 | 10 | rowsb = p; | |
114 | 10 | colsb = ds->extrarow? m+1: m; | |
115 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
10 | if (format == PETSC_VIEWER_ASCII_MATLAB) { |
116 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer,"%% Size = %" PetscInt_FMT " %" PetscInt_FMT "\n",rowsa,colsa)); | |
117 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer,"zzz = zeros(%" PetscInt_FMT ",3);\n",2*ds->n)); | |
118 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer,"zzz = [\n")); | |
119 | ✗ | for (i=0;i<PetscMin(rowsa,colsa);i++) PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT " %18.16e\n",i+1,i+1,(double)T[i])); | |
120 | ✗ | for (i=0;i<k;i++) PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT " %18.16e\n",i+1,k+1,(double)T[i+ds->ld])); | |
121 | ✗ | if (n>m) { /* A lower bidiagonal */ | |
122 | ✗ | for (i=k;i<rowsa-1;i++) PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT " %18.16e\n",i+2,i+1,(double)T[i+ds->ld])); | |
123 | } else { /* A (square) upper bidiagonal */ | ||
124 | ✗ | for (i=k;i<colsa-1;i++) PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT " %18.16e\n",i+1,i+2,(double)T[i+ds->ld])); | |
125 | } | ||
126 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_T])); | |
127 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer,"%% Size = %" PetscInt_FMT " %" PetscInt_FMT "\n",rowsb,colsb)); | |
128 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer,"zzz = zeros(%" PetscInt_FMT ",3);\n",2*ds->n)); | |
129 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer,"zzz = [\n")); | |
130 | ✗ | for (i=0;i<rowsb;i++) PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT " %18.16e\n",i+1,i+1,(double)S[i])); | |
131 | ✗ | for (i=0;i<colsb-1;i++) { | |
132 | ✗ | r = PetscMax(i+2,ds->k+1); | |
133 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT " %18.16e\n",i+1,r,(double)T[i+2*ds->ld])); | |
134 | } | ||
135 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_D])); | |
136 | } else { | ||
137 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscViewerASCIIPrintf(viewer,"Matrix %s =\n",DSMatName[DS_MAT_T])); |
138 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
110 | for (i=0;i<rowsa;i++) { |
139 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1100 | for (j=0;j<colsa;j++) { |
140 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1000 | if (i==j) value = T[i]; |
141 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
900 | else if (i<ds->k && j==ds->k) value = T[i+ds->ld]; |
142 |
1/6✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
900 | else if (n>m && i==j+1 && i>ds->k) value = T[j+ds->ld]; |
143 |
4/6✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 10 times.
✗ Branch 5 not taken.
|
900 | else if (n<=m && i+1==j && i>=ds->k) value = T[i+ds->ld]; |
144 | else value = 0.0; | ||
145 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1000 | PetscCall(PetscViewerASCIIPrintf(viewer," %18.16e ",(double)value)); |
146 | } | ||
147 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
100 | PetscCall(PetscViewerASCIIPrintf(viewer,"\n")); |
148 | } | ||
149 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscViewerASCIIPrintf(viewer,"Matrix %s =\n",DSMatName[DS_MAT_D])); |
150 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
110 | for (i=0;i<rowsb;i++) { |
151 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1100 | for (j=0;j<colsb;j++) { |
152 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1000 | if (i==j) value = S[i]; |
153 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
900 | else if (i<ds->k && j==ds->k) value = T[PetscMin(i,j)+2*ds->ld]; |
154 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
|
900 | else if (i+1==j && i>=ds->k) value = T[i+2*ds->ld]; |
155 | else value = 0.0; | ||
156 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1000 | PetscCall(PetscViewerASCIIPrintf(viewer," %18.16e ",(double)value)); |
157 | } | ||
158 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
100 | PetscCall(PetscViewerASCIIPrintf(viewer,"\n")); |
159 | } | ||
160 | } | ||
161 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); |
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.
|
10 | PetscCall(PetscViewerFlush(viewer)); |
163 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T)); |
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.
|
10 | PetscCall(DSRestoreArrayReal(ds,DS_MAT_D,&S)); |
165 | } else { | ||
166 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(DSViewMat(ds,viewer,DS_MAT_A)); |
167 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(DSViewMat(ds,viewer,DS_MAT_B)); |
168 | } | ||
169 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
130 | if (ds->state>DS_STATE_INTERMEDIATE) { |
170 | ✗ | PetscCall(DSViewMat(ds,viewer,DS_MAT_X)); | |
171 | ✗ | PetscCall(DSViewMat(ds,viewer,DS_MAT_U)); | |
172 | ✗ | PetscCall(DSViewMat(ds,viewer,DS_MAT_V)); | |
173 | } | ||
174 |
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.
|
26 | PetscFunctionReturn(PETSC_SUCCESS); |
175 | } | ||
176 | |||
177 | 170 | static PetscErrorCode DSVectors_GSVD(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm) | |
178 | { | ||
179 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
170 | PetscFunctionBegin; |
180 |
1/3✗ Branch 0 not taken.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
|
170 | switch (mat) { |
181 | ✗ | case DS_MAT_U: | |
182 | case DS_MAT_V: | ||
183 | ✗ | if (rnorm) *rnorm = 0.0; | |
184 | break; | ||
185 | case DS_MAT_X: | ||
186 | break; | ||
187 | ✗ | default: | |
188 | ✗ | SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter"); | |
189 | } | ||
190 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
34 | PetscFunctionReturn(PETSC_SUCCESS); |
191 | } | ||
192 | |||
193 | 6331 | static PetscErrorCode DSSort_GSVD(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k) | |
194 | { | ||
195 | 6331 | DS_GSVD *ctx = (DS_GSVD*)ds->data; | |
196 | 6331 | PetscInt t,l,ld=ds->ld,i,*perm,*perm2; | |
197 | 6331 | PetscReal *T=NULL,*D=NULL,*eig; | |
198 | 6331 | PetscScalar *A=NULL,*B=NULL; | |
199 | 6331 | PetscBool compact=ds->compact; | |
200 | |||
201 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
6331 | PetscFunctionBegin; |
202 |
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.
|
6331 | if (!ds->sc) PetscFunctionReturn(PETSC_SUCCESS); |
203 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6331 | PetscCheck(ctx->m,PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the other dimensions with DSGSVDSetDimensions()"); |
204 | 6331 | l = ds->l; | |
205 | 6331 | t = ds->t; | |
206 | 6331 | perm = ds->perm; | |
207 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(PetscMalloc2(t,&eig,t,&perm2)); |
208 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
6331 | if (compact) { |
209 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6141 | PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T)); |
210 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6141 | PetscCall(DSGetArrayReal(ds,DS_MAT_D,&D)); |
211 |
3/4✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
63763 | for (i=0;i<t;i++) eig[i] = (D[i]==0)?PETSC_INFINITY:T[i]/D[i]; |
212 | } else { | ||
213 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
190 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_A],&A)); |
214 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
190 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_B],&B)); |
215 |
3/4✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
2300 | for (i=0;i<t;i++) eig[i] = (B[i+i*ld]==0)?PETSC_INFINITY:PetscRealPart(A[i+i*ld])/PetscRealPart(B[i*(1+ld)]); |
216 | } | ||
217 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(DSSortEigenvaluesReal_Private(ds,eig,perm)); |
218 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(PetscArraycpy(perm2,perm,t)); |
219 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
62688 | for (i=l;i<t;i++) wr[i] = eig[perm[i]]; |
220 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
6331 | if (compact) { |
221 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6141 | PetscCall(PetscArraycpy(eig,T,t)); |
222 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
60388 | for (i=l;i<t;i++) T[i] = eig[perm[i]]; |
223 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6141 | PetscCall(PetscArraycpy(eig,D,t)); |
224 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
60388 | for (i=l;i<t;i++) D[i] = eig[perm[i]]; |
225 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6141 | PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T)); |
226 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6141 | PetscCall(DSRestoreArrayReal(ds,DS_MAT_D,&D)); |
227 | } else { | ||
228 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2300 | for (i=l;i<t;i++) eig[i] = PetscRealPart(A[i*(1+ld)]); |
229 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2300 | for (i=l;i<t;i++) A[i*(1+ld)] = eig[perm[i]]; |
230 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2300 | for (i=l;i<t;i++) eig[i] = PetscRealPart(B[i*(1+ld)]); |
231 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2300 | for (i=l;i<t;i++) B[i*(1+ld)] = eig[perm[i]]; |
232 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
190 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_A],&A)); |
233 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
190 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_B],&B)); |
234 | } | ||
235 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(DSPermuteColumns_Private(ds,l,t,ds->n,DS_MAT_U,perm2)); |
236 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(PetscArraycpy(perm2,perm,t)); |
237 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(DSPermuteColumns_Private(ds,l,t,ctx->m,DS_MAT_X,perm2)); |
238 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(DSPermuteColumns_Private(ds,l,t,ctx->p,DS_MAT_V,perm)); |
239 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(PetscFree2(eig,perm2)); |
240 |
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.
|
1254 | PetscFunctionReturn(PETSC_SUCCESS); |
241 | } | ||
242 | |||
243 | 6111 | static PetscErrorCode DSUpdateExtraRow_GSVD(DS ds) | |
244 | { | ||
245 | 6111 | DS_GSVD *ctx = (DS_GSVD*)ds->data; | |
246 | 6111 | PetscInt i; | |
247 | 6111 | PetscBLASInt n=0,m=0,ld=0; | |
248 | 6111 | const PetscScalar *U,*V; | |
249 | 6111 | PetscReal *T,*e,*f,alpha,beta,betah; | |
250 | |||
251 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
6111 | PetscFunctionBegin; |
252 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6111 | PetscCheck(ctx->m,PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the other dimensions with DSGSVDSetDimensions()"); |
253 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6111 | PetscCheck(ds->compact,PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented for non-compact storage"); |
254 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6111 | PetscCall(PetscBLASIntCast(ds->n,&n)); |
255 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6111 | PetscCall(PetscBLASIntCast(ctx->m,&m)); |
256 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6111 | PetscCall(PetscBLASIntCast(ds->ld,&ld)); |
257 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6111 | PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T)); |
258 | 6111 | e = T+ld; | |
259 | 6111 | f = T+2*ld; | |
260 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6111 | PetscCall(MatDenseGetArrayRead(ds->omat[DS_MAT_U],&U)); |
261 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6111 | PetscCall(MatDenseGetArrayRead(ds->omat[DS_MAT_V],&V)); |
262 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
6111 | if (n<=m) { /* upper variant, A is square upper bidiagonal */ |
263 | 2488 | beta = e[m-1]; /* in compact, we assume all entries are zero except the last one */ | |
264 | 2488 | betah = f[m-1]; | |
265 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
25084 | for (i=0;i<m;i++) { |
266 | 22596 | e[i] = PetscRealPart(beta*U[m-1+i*ld]); | |
267 | 22596 | f[i] = PetscRealPart(betah*V[m-1+i*ld]); | |
268 | } | ||
269 | } else { /* lower variant, A is (m+1)xm lower bidiagonal */ | ||
270 | 3623 | alpha = T[m]; | |
271 | 3623 | betah = f[m-1]; | |
272 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
38349 | for (i=0;i<m;i++) { |
273 | 34726 | e[i] = PetscRealPart(alpha*U[m+i*ld]); | |
274 | 34726 | f[i] = PetscRealPart(betah*V[m-1+i*ld]); | |
275 | } | ||
276 | 3623 | T[m] = PetscRealPart(alpha*U[m+m*ld]); | |
277 | } | ||
278 | 6111 | ds->k = m; | |
279 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6111 | PetscCall(MatDenseRestoreArrayRead(ds->omat[DS_MAT_U],&U)); |
280 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6111 | PetscCall(MatDenseRestoreArrayRead(ds->omat[DS_MAT_V],&V)); |
281 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6111 | PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T)); |
282 |
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.
|
1210 | PetscFunctionReturn(PETSC_SUCCESS); |
283 | } | ||
284 | |||
285 | 5986 | static PetscErrorCode DSTruncate_GSVD(DS ds,PetscInt n,PetscBool trim) | |
286 | { | ||
287 | 5986 | DS_GSVD *ctx = (DS_GSVD*)ds->data; | |
288 | 5986 | PetscScalar *U; | |
289 | 5986 | PetscReal *T; | |
290 | 5986 | PetscInt i,m=ctx->m,ld=ds->ld; | |
291 | 5986 | PetscBool lower=(ds->n>ctx->m)?PETSC_TRUE:PETSC_FALSE; | |
292 | |||
293 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
5986 | PetscFunctionBegin; |
294 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
5986 | PetscCheck(ds->compact,PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented for non-compact storage"); |
295 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
5986 | if (trim) { |
296 | 508 | ds->l = 0; | |
297 | 508 | ds->k = 0; | |
298 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
508 | ds->n = lower? n+1: n; |
299 | 508 | ctx->m = n; | |
300 | 508 | ctx->p = n; | |
301 | 508 | ds->t = ds->n; /* truncated length equal to the new dimension */ | |
302 | 508 | ctx->tm = ctx->m; /* must also keep the previous dimension of X */ | |
303 | 508 | ctx->tp = ctx->p; /* must also keep the previous dimension of V */ | |
304 | } else { | ||
305 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
5478 | if (lower) { |
306 | /* move value of diagonal element of arrow (alpha) */ | ||
307 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3193 | PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T)); |
308 | 3193 | T[n] = T[m]; | |
309 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3193 | PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T)); |
310 | /* copy last column of U so that it updates the next initial vector of U1 */ | ||
311 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3193 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_U],&U)); |
312 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
36764 | for (i=0;i<=m;i++) U[i+n*ld] = U[i+m*ld]; |
313 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3193 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_U],&U)); |
314 | } | ||
315 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
5478 | ds->k = (ds->extrarow)? n: 0; |
316 | 5478 | ds->t = ds->n; /* truncated length equal to previous dimension */ | |
317 | 5478 | ctx->tm = ctx->m; /* must also keep the previous dimension of X */ | |
318 | 5478 | ctx->tp = ctx->p; /* must also keep the previous dimension of V */ | |
319 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
5478 | ds->n = lower? n+1: n; |
320 | 5478 | ctx->m = n; | |
321 | 5478 | ctx->p = n; | |
322 | } | ||
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.
|
1185 | PetscFunctionReturn(PETSC_SUCCESS); |
324 | } | ||
325 | |||
326 | 6261 | static PetscErrorCode DSSwitchFormat_GSVD(DS ds) | |
327 | { | ||
328 | 6261 | DS_GSVD *ctx = (DS_GSVD*)ds->data; | |
329 | 6261 | PetscReal *T,*D; | |
330 | 6261 | PetscScalar *A,*B; | |
331 | 6261 | PetscInt i,n=ds->n,k=ds->k,ld=ds->ld,m=ctx->m; | |
332 | |||
333 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
6261 | PetscFunctionBegin; |
334 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6261 | PetscCheck(ctx->m,PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the other dimensions with DSGSVDSetDimensions()"); |
335 | /* switch from compact (arrow) to dense storage */ | ||
336 | /* bidiagonal associated to B is stored in D and T+2*ld */ | ||
337 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6261 | PetscCall(MatDenseGetArrayWrite(ds->omat[DS_MAT_A],&A)); |
338 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6261 | PetscCall(MatDenseGetArrayWrite(ds->omat[DS_MAT_B],&B)); |
339 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6261 | PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T)); |
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.
|
6261 | PetscCall(DSGetArrayReal(ds,DS_MAT_D,&D)); |
341 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6261 | PetscCall(PetscArrayzero(A,ld*ld)); |
342 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6261 | PetscCall(PetscArrayzero(B,ld*ld)); |
343 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
30939 | for (i=0;i<k;i++) { |
344 | 24678 | A[i+i*ld] = T[i]; | |
345 | 24678 | A[i+k*ld] = T[i+ld]; | |
346 | 24678 | B[i+i*ld] = D[i]; | |
347 | 24678 | B[i+k*ld] = T[i+2*ld]; | |
348 | } | ||
349 | /* B is upper bidiagonal */ | ||
350 | 6261 | B[k+k*ld] = D[k]; | |
351 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
33914 | for (i=k+1;i<m;i++) { |
352 | 27653 | B[i+i*ld] = D[i]; | |
353 | 27653 | B[i-1+i*ld] = T[i-1+2*ld]; | |
354 | } | ||
355 | /* A can be upper (square) or lower bidiagonal */ | ||
356 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
40105 | for (i=k;i<m;i++) A[i+i*ld] = T[i]; |
357 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
26719 | if (n>m) for (i=k;i<m;i++) A[i+1+i*ld] = T[i+ld]; |
358 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
13386 | else for (i=k+1;i<m;i++) A[i-1+i*ld] = T[i-1+ld]; |
359 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6261 | PetscCall(MatDenseRestoreArrayWrite(ds->omat[DS_MAT_A],&A)); |
360 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6261 | PetscCall(MatDenseRestoreArrayWrite(ds->omat[DS_MAT_B],&B)); |
361 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6261 | PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T)); |
362 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6261 | PetscCall(DSRestoreArrayReal(ds,DS_MAT_D,&D)); |
363 |
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.
|
1240 | PetscFunctionReturn(PETSC_SUCCESS); |
364 | } | ||
365 | |||
366 | /* | ||
367 | Compact format is used when [A;B] has orthonormal columns. | ||
368 | In this case R=I and the GSVD of (A,B) is the CS decomposition | ||
369 | */ | ||
370 | 6331 | static PetscErrorCode DSSolve_GSVD(DS ds,PetscScalar *wr,PetscScalar *wi) | |
371 | { | ||
372 | 6331 | DS_GSVD *ctx = (DS_GSVD*)ds->data; | |
373 | 6331 | PetscInt i,j; | |
374 | 6331 | PetscBLASInt n1,m1,info,lc = 0,n = 0,m = 0,p = 0,p1,l,k,q,ld,off,lwork,r; | |
375 | 6331 | PetscScalar *A,*B,*X,*U,*V,sone=1.0,smone=-1.0; | |
376 | 6331 | PetscReal *alpha,*beta,*T,*D; | |
377 | #if !defined(SLEPC_MISSING_LAPACK_GGSVD3) | ||
378 | 6331 | PetscScalar a,dummy; | |
379 | 6331 | PetscReal rdummy; | |
380 | 6331 | PetscBLASInt idummy; | |
381 | #endif | ||
382 | |||
383 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
6331 | PetscFunctionBegin; |
384 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6331 | PetscCheck(ctx->m,PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the other dimensions with DSGSVDSetDimensions()"); |
385 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(PetscBLASIntCast(ds->n,&m)); |
386 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(PetscBLASIntCast(ctx->m,&n)); |
387 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(PetscBLASIntCast(ctx->p,&p)); |
388 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(PetscBLASIntCast(ds->l,&lc)); |
389 |
3/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
6331 | PetscCheck(ds->compact || lc==0,PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DSGSVD with non-compact format does not support locking"); |
390 | /* In compact storage B is always nxn and A can be either nxn or (n+1)xn */ | ||
391 |
6/10✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 10 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 10 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
6331 | PetscCheck(!ds->compact || (p==n && (m==p || m==p+1)),PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Dimensions not supported in compact format"); |
392 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(PetscBLASIntCast(ds->ld,&ld)); |
393 | 6331 | n1 = n-lc; /* n1 = size of leading block, excl. locked + size of trailing block */ | |
394 | 6331 | m1 = m-lc; | |
395 | 6331 | p1 = p-lc; | |
396 | 6331 | off = lc+lc*ld; | |
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.
|
6331 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_A],&A)); |
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.
|
6331 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_B],&B)); |
399 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_X],&X)); |
400 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_U],&U)); |
401 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_V],&V)); |
402 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(PetscArrayzero(X,ld*ld)); |
403 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
9706 | for (i=0;i<lc;i++) X[i+i*ld] = 1.0; |
404 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(PetscArrayzero(U,ld*ld)); |
405 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
9706 | for (i=0;i<lc;i++) U[i+i*ld] = 1.0; |
406 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(PetscArrayzero(V,ld*ld)); |
407 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
9706 | for (i=0;i<lc;i++) V[i+i*ld] = 1.0; |
408 |
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.
|
6331 | if (ds->compact) PetscCall(DSSwitchFormat_GSVD(ds)); |
409 | |||
410 | #if !defined(SLEPC_MISSING_LAPACK_GGSVD3) | ||
411 | /* workspace query and memory allocation */ | ||
412 | 6331 | lwork = -1; | |
413 | #if !defined (PETSC_USE_COMPLEX) | ||
414 |
10/20✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
✓ Branch 12 taken 1 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 1 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 1 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
3521 | PetscCallBLAS("LAPACKggsvd3",LAPACKggsvd3_("U","V","Q",&m1,&n1,&p1,&k,&l,&dummy,&ld,&dummy,&ld,&rdummy,&rdummy,&dummy,&ld,&dummy,&ld,&dummy,&ld,&a,&lwork,&idummy,&info)); |
415 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
3521 | PetscCall(PetscBLASIntCast((PetscInt)a,&lwork)); |
416 | #else | ||
417 |
10/20✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
✓ Branch 12 taken 1 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 1 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 1 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
2810 | PetscCallBLAS("LAPACKggsvd3",LAPACKggsvd3_("U","V","Q",&m1,&n1,&p1,&k,&l,&dummy,&ld,&dummy,&ld,&rdummy,&rdummy,&dummy,&ld,&dummy,&ld,&dummy,&ld,&a,&lwork,&rdummy,&idummy,&info)); |
418 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
2810 | PetscCall(PetscBLASIntCast((PetscInt)PetscRealPart(a),&lwork)); |
419 | #endif | ||
420 | |||
421 | #if !defined (PETSC_USE_COMPLEX) | ||
422 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
3521 | PetscCall(DSAllocateWork_Private(ds,lwork,2*ds->ld,ds->ld)); |
423 | 3521 | alpha = ds->rwork; | |
424 | 3521 | beta = ds->rwork+ds->ld; | |
425 |
10/20✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
✓ Branch 12 taken 1 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 1 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 1 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
3521 | PetscCallBLAS("LAPACKggsvd3",LAPACKggsvd3_("U","V","Q",&m1,&n1,&p1,&k,&l,A+off,&ld,B+off,&ld,alpha,beta,U+off,&ld,V+off,&ld,X+off,&ld,ds->work,&lwork,ds->iwork,&info)); |
426 | #else | ||
427 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
2810 | PetscCall(DSAllocateWork_Private(ds,lwork,4*ds->ld,ds->ld)); |
428 | 2810 | alpha = ds->rwork+2*ds->ld; | |
429 | 2810 | beta = ds->rwork+3*ds->ld; | |
430 |
10/20✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
✓ Branch 12 taken 1 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 1 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 1 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
2810 | PetscCallBLAS("LAPACKggsvd3",LAPACKggsvd3_("U","V","Q",&m1,&n1,&p1,&k,&l,A+off,&ld,B+off,&ld,alpha,beta,U+off,&ld,V+off,&ld,X+off,&ld,ds->work,&lwork,ds->rwork,ds->iwork,&info)); |
431 | #endif | ||
432 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6331 | SlepcCheckLapackInfo("ggsvd3",info); |
433 | |||
434 | #else /* defined(SLEPC_MISSING_LAPACK_GGSVD3) */ | ||
435 | |||
436 | lwork = PetscMax(PetscMax(3*n,m),p)+n; | ||
437 | #if !defined (PETSC_USE_COMPLEX) | ||
438 | PetscCall(DSAllocateWork_Private(ds,lwork,2*ds->ld,ds->ld)); | ||
439 | alpha = ds->rwork; | ||
440 | beta = ds->rwork+ds->ld; | ||
441 | PetscCallBLAS("LAPACKggsvd",LAPACKggsvd_("U","V","Q",&m1,&n1,&p1,&k,&l,A+off,&ld,B+off,&ld,alpha,beta,U+off,&ld,V+off,&ld,X+off,&ld,ds->work,ds->iwork,&info)); | ||
442 | #else | ||
443 | PetscCall(DSAllocateWork_Private(ds,lwork,4*ds->ld,ds->ld)); | ||
444 | alpha = ds->rwork+2*ds->ld; | ||
445 | beta = ds->rwork+3*ds->ld; | ||
446 | PetscCallBLAS("LAPACKggsvd",LAPACKggsvd_("U","V","Q",&m1,&n1,&p1,&k,&l,A+off,&ld,B+off,&ld,alpha,beta,U+off,&ld,V+off,&ld,X+off,&ld,ds->work,ds->rwork,ds->iwork,&info)); | ||
447 | #endif | ||
448 | SlepcCheckLapackInfo("ggsvd",info); | ||
449 | |||
450 | #endif | ||
451 | |||
452 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6331 | PetscCheck(k+l>=n1,PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"The rank deficient case not supported yet"); |
453 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
6331 | if (ds->compact) { |
454 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6141 | PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T)); |
455 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6141 | PetscCall(DSGetArrayReal(ds,DS_MAT_D,&D)); |
456 | /* R is the identity matrix (except the sign) */ | ||
457 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
60388 | for (i=lc;i<n;i++) { |
458 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
54247 | if (PetscRealPart(A[i+i*ld])<0.0) { /* scale column i */ |
459 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
320242 | for (j=lc;j<n;j++) X[j+i*ld] = -X[j+i*ld]; |
460 | } | ||
461 | } | ||
462 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6141 | PetscCall(PetscArrayzero(T+ld,m-1)); |
463 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6141 | PetscCall(PetscArrayzero(T+2*ld,n-1)); |
464 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
60388 | for (i=lc;i<n;i++) { |
465 | 54247 | T[i] = alpha[i-lc]; | |
466 | 54247 | D[i] = beta[i-lc]; | |
467 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
54247 | if (D[i]==0.0) wr[i] = PETSC_INFINITY; |
468 | 54247 | else wr[i] = T[i]/D[i]; | |
469 | } | ||
470 | 6141 | ds->t = n; | |
471 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6141 | PetscCall(DSRestoreArrayReal(ds,DS_MAT_D,&D)); |
472 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6141 | PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T)); |
473 | } else { | ||
474 | /* X = X*inv(R) */ | ||
475 | 190 | q = PetscMin(m,n); | |
476 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
190 | PetscCallBLAS("BLAStrsm",BLAStrsm_("R","U","N","N",&n,&q,&sone,A,&ld,X,&ld)); |
477 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
190 | if (m<n) { |
478 | 40 | r = n-m; | |
479 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
40 | PetscCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&r,&m,&sone,X,&ld,A,&ld,&smone,X+m*ld,&ld)); |
480 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
40 | PetscCallBLAS("BLAStrsm",BLAStrsm_("R","U","N","N",&n,&r,&sone,B+m*ld,&ld,X+m*ld,&ld)); |
481 | } | ||
482 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
190 | if (k>0) { |
483 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
130 | for (i=k;i<PetscMin(m,k+l);i++) { |
484 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
110 | PetscCall(PetscArraycpy(X+(i-k)*ld,X+i*ld,ld)); |
485 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
110 | PetscCall(PetscArraycpy(U+(i-k)*ld,U+i*ld,ld)); |
486 | } | ||
487 | } | ||
488 | /* singular values */ | ||
489 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
190 | PetscCall(PetscArrayzero(A,ld*ld)); |
490 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
190 | PetscCall(PetscArrayzero(B,ld*ld)); |
491 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2300 | for (j=k;j<PetscMin(m,k+l);j++) { |
492 | 2110 | A[(j-k)*(1+ld)] = alpha[j]; | |
493 | 2110 | B[(j-k)*(1+ld)] = beta[j]; | |
494 | 2110 | wr[j-k] = alpha[j]/beta[j]; | |
495 | } | ||
496 | 190 | ds->t = PetscMin(m,k+l)-k; /* set number of computed values */ | |
497 | } | ||
498 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_A],&A)); |
499 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_B],&B)); |
500 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_X],&X)); |
501 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_U],&U)); |
502 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6331 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_V],&V)); |
503 |
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.
|
1254 | PetscFunctionReturn(PETSC_SUCCESS); |
504 | } | ||
505 | |||
506 | 240 | static PetscErrorCode DSCond_GSVD(DS ds,PetscReal *cond) | |
507 | { | ||
508 | 240 | DS_GSVD *ctx = (DS_GSVD*)ds->data; | |
509 | 240 | PetscBLASInt lwork,lrwork=0,info,m,n,p,ld; | |
510 | 240 | PetscScalar *A,*work; | |
511 | 240 | const PetscScalar *M; | |
512 | 240 | PetscReal *sigma,conda,condb; | |
513 | #if defined(PETSC_USE_COMPLEX) | ||
514 | 120 | PetscReal *rwork; | |
515 | #endif | ||
516 | |||
517 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
240 | PetscFunctionBegin; |
518 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
240 | PetscCall(PetscBLASIntCast(ds->n,&m)); |
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.
|
240 | PetscCall(PetscBLASIntCast(ctx->m,&n)); |
520 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
240 | PetscCall(PetscBLASIntCast(ctx->p,&p)); |
521 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
240 | PetscCall(PetscBLASIntCast(ds->ld,&ld)); |
522 | 240 | lwork = 5*n; | |
523 | #if defined(PETSC_USE_COMPLEX) | ||
524 | 120 | lrwork = 5*n; | |
525 | #endif | ||
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.
|
240 | PetscCall(DSAllocateWork_Private(ds,ld*n+lwork,n+lrwork,0)); |
527 | 240 | A = ds->work; | |
528 | 240 | work = ds->work+ld*n; | |
529 | 240 | sigma = ds->rwork; | |
530 | #if defined(PETSC_USE_COMPLEX) | ||
531 | 120 | rwork = ds->rwork+n; | |
532 | #endif | ||
533 |
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.
|
240 | if (ds->compact) PetscCall(DSSwitchFormat_GSVD(ds)); |
534 | |||
535 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
240 | PetscCall(MatDenseGetArrayRead(ds->omat[DS_MAT_A],&M)); |
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.
|
240 | PetscCall(PetscArraycpy(A,M,ld*n)); |
537 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
240 | PetscCall(MatDenseRestoreArrayRead(ds->omat[DS_MAT_A],&M)); |
538 | #if defined(PETSC_USE_COMPLEX) | ||
539 |
10/20✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
✓ Branch 12 taken 1 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 1 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 1 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
120 | PetscCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&m,&n,A,&ld,sigma,NULL,&ld,NULL,&ld,work,&lwork,rwork,&info)); |
540 | #else | ||
541 |
10/20✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
✓ Branch 12 taken 1 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 1 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 1 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
120 | PetscCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&m,&n,A,&ld,sigma,NULL,&ld,NULL,&ld,work,&lwork,&info)); |
542 | #endif | ||
543 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
240 | SlepcCheckLapackInfo("gesvd",info); |
544 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
240 | conda = sigma[0]/sigma[PetscMin(m,n)-1]; |
545 | |||
546 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
240 | PetscCall(MatDenseGetArrayRead(ds->omat[DS_MAT_B],&M)); |
547 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
240 | PetscCall(PetscArraycpy(A,M,ld*n)); |
548 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
240 | PetscCall(MatDenseRestoreArrayRead(ds->omat[DS_MAT_B],&M)); |
549 | #if defined(PETSC_USE_COMPLEX) | ||
550 |
10/20✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
✓ Branch 12 taken 1 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 1 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 1 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
120 | PetscCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&p,&n,A,&ld,sigma,NULL,&ld,NULL,&ld,work,&lwork,rwork,&info)); |
551 | #else | ||
552 |
10/20✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
✓ Branch 12 taken 1 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 1 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 1 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
120 | PetscCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&p,&n,A,&ld,sigma,NULL,&ld,NULL,&ld,work,&lwork,&info)); |
553 | #endif | ||
554 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
240 | SlepcCheckLapackInfo("gesvd",info); |
555 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
240 | condb = sigma[0]/sigma[PetscMin(p,n)-1]; |
556 | |||
557 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
240 | *cond = PetscMax(conda,condb); |
558 |
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.
|
240 | PetscFunctionReturn(PETSC_SUCCESS); |
559 | } | ||
560 | |||
561 | #if !defined(PETSC_HAVE_MPIUNI) | ||
562 | 300 | static PetscErrorCode DSSynchronize_GSVD(DS ds,PetscScalar eigr[],PetscScalar eigi[]) | |
563 | { | ||
564 | 300 | DS_GSVD *ctx = (DS_GSVD*)ds->data; | |
565 | 300 | PetscInt ld=ds->ld,l=ds->l,k=0,kr=0; | |
566 | 300 | PetscMPIInt m,rank,off=0,size,n,ldn,ld3; | |
567 | 300 | PetscScalar *A,*U,*V,*X; | |
568 | 300 | PetscReal *T; | |
569 | |||
570 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
300 | PetscFunctionBegin; |
571 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
300 | PetscCall(PetscMPIIntCast(ctx->m,&m)); |
572 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
300 | if (ds->compact) kr = 3*ld; |
573 | 50 | else k = 2*(m-l)*ld; | |
574 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
300 | if (ds->state>DS_STATE_RAW) k += 3*(m-l)*ld; |
575 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
300 | if (eigr) k += m-l; |
576 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
300 | PetscCall(DSAllocateWork_Private(ds,k+kr,0,0)); |
577 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
300 | PetscCall(PetscMPIIntCast(k*sizeof(PetscScalar)+kr*sizeof(PetscReal),&size)); |
578 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
300 | PetscCall(PetscMPIIntCast(m-l,&n)); |
579 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
300 | PetscCall(PetscMPIIntCast(ld*(m-l),&ldn)); |
580 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
300 | PetscCall(PetscMPIIntCast(3*ld,&ld3)); |
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.
|
300 | if (ds->compact) PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T)); |
582 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | else PetscCall(MatDenseGetArray(ds->omat[DS_MAT_A],&A)); |
583 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
300 | if (ds->state>DS_STATE_RAW) { |
584 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
300 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_U],&U)); |
585 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
300 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_V],&V)); |
586 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
300 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_X],&X)); |
587 | } | ||
588 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
300 | PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank)); |
589 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
300 | if (!rank) { |
590 |
16/30✓ 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 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.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
145 | if (ds->compact) PetscCallMPI(MPI_Pack(T,ld3,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
591 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
20 | else PetscCallMPI(MPI_Pack(A+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
592 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
145 | if (ds->state>DS_STATE_RAW) { |
593 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
145 | PetscCallMPI(MPI_Pack(U+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
594 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
145 | PetscCallMPI(MPI_Pack(V+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
595 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
145 | PetscCallMPI(MPI_Pack(X+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
596 | } | ||
597 |
15/30✓ 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 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.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
145 | if (eigr) PetscCallMPI(MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
598 | } | ||
599 |
15/30✓ 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 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.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
600 | PetscCallMPI(MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds))); |
600 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
300 | if (rank) { |
601 |
16/30✓ 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 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.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
155 | if (ds->compact) PetscCallMPI(MPI_Unpack(ds->work,size,&off,T,ld3,MPIU_REAL,PetscObjectComm((PetscObject)ds))); |
602 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
30 | else PetscCallMPI(MPI_Unpack(ds->work,size,&off,A+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds))); |
603 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
155 | if (ds->state>DS_STATE_RAW) { |
604 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
155 | PetscCallMPI(MPI_Unpack(ds->work,size,&off,U+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds))); |
605 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
155 | PetscCallMPI(MPI_Unpack(ds->work,size,&off,V+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds))); |
606 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
155 | PetscCallMPI(MPI_Unpack(ds->work,size,&off,X+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds))); |
607 | } | ||
608 |
15/30✓ 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 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.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
155 | if (eigr) PetscCallMPI(MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds))); |
609 | } | ||
610 |
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.
|
300 | if (ds->compact) PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T)); |
611 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | else PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_A],&A)); |
612 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
300 | if (ds->state>DS_STATE_RAW) { |
613 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
300 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_U],&U)); |
614 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
300 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_V],&V)); |
615 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
300 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_X],&X)); |
616 | } | ||
617 |
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.
|
60 | PetscFunctionReturn(PETSC_SUCCESS); |
618 | } | ||
619 | #endif | ||
620 | |||
621 | 18848 | static PetscErrorCode DSMatGetSize_GSVD(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols) | |
622 | { | ||
623 | 18848 | DS_GSVD *ctx = (DS_GSVD*)ds->data; | |
624 | |||
625 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
18848 | PetscFunctionBegin; |
626 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
18848 | PetscCheck(ctx->m,PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the other dimensions with DSGSVDSetDimensions()"); |
627 |
5/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 10 times.
✓ Branch 6 taken 10 times.
✗ Branch 7 not taken.
|
18848 | switch (t) { |
628 | 190 | case DS_MAT_A: | |
629 | 190 | *rows = ds->n; | |
630 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
190 | *cols = ds->extrarow? ctx->m+1: ctx->m; |
631 | 190 | break; | |
632 | 190 | case DS_MAT_B: | |
633 | 190 | *rows = ctx->p; | |
634 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
190 | *cols = ds->extrarow? ctx->m+1: ctx->m; |
635 | 190 | break; | |
636 | ✗ | case DS_MAT_T: | |
637 | ✗ | *rows = ds->n; | |
638 | ✗ | *cols = PetscDefined(USE_COMPLEX)? 2: 3; | |
639 | ✗ | break; | |
640 | ✗ | case DS_MAT_D: | |
641 | ✗ | *rows = ctx->p; | |
642 | ✗ | *cols = 1; | |
643 | ✗ | break; | |
644 | 6156 | case DS_MAT_U: | |
645 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
6156 | *rows = ds->state==DS_STATE_TRUNCATED? ds->t: ds->n; |
646 | 6156 | *cols = ds->n; | |
647 | 6156 | break; | |
648 | 6156 | case DS_MAT_V: | |
649 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
6156 | *rows = ds->state==DS_STATE_TRUNCATED? ctx->tp: ctx->p; |
650 | 6156 | *cols = ctx->p; | |
651 | 6156 | break; | |
652 | 6156 | case DS_MAT_X: | |
653 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
6156 | *rows = ds->state==DS_STATE_TRUNCATED? ctx->tm: ctx->m; |
654 | 6156 | *cols = ctx->m; | |
655 | 6156 | break; | |
656 | ✗ | default: | |
657 | ✗ | SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid t parameter"); | |
658 | } | ||
659 |
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.
|
3733 | PetscFunctionReturn(PETSC_SUCCESS); |
660 | } | ||
661 | |||
662 | 6331 | static PetscErrorCode DSGSVDSetDimensions_GSVD(DS ds,PetscInt m,PetscInt p) | |
663 | { | ||
664 | 6331 | DS_GSVD *ctx = (DS_GSVD*)ds->data; | |
665 | |||
666 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
6331 | PetscFunctionBegin; |
667 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6331 | DSCheckAlloc(ds,1); |
668 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
6331 | if (m == PETSC_DETERMINE) { |
669 | ✗ | ctx->m = ds->ld; | |
670 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
6331 | } else if (m != PETSC_CURRENT) { |
671 |
2/6✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
6331 | PetscCheck(m>0 && m<=ds->ld,PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of m. Must be between 1 and ld"); |
672 | 6331 | ctx->m = m; | |
673 | } | ||
674 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
6331 | if (p == PETSC_DETERMINE) { |
675 | 30 | ctx->p = ds->n; | |
676 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
6301 | } else if (p != PETSC_CURRENT) { |
677 |
2/6✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
6301 | PetscCheck(p>0 && p<=ds->ld,PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of p. Must be between 1 and ld"); |
678 | 6301 | ctx->p = p; | |
679 | } | ||
680 |
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.
|
1254 | PetscFunctionReturn(PETSC_SUCCESS); |
681 | } | ||
682 | |||
683 | /*@ | ||
684 | DSGSVDSetDimensions - Sets the number of columns and rows for a DSGSVD. | ||
685 | |||
686 | Logically Collective | ||
687 | |||
688 | Input Parameters: | ||
689 | + ds - the direct solver context | ||
690 | . m - the number of columns | ||
691 | - p - the number of rows for the second matrix (B) | ||
692 | |||
693 | Notes: | ||
694 | This call is complementary to DSSetDimensions(), to provide two dimensions | ||
695 | that are specific to this DS type. The number of rows for the first matrix (A) | ||
696 | is set by DSSetDimensions(). | ||
697 | |||
698 | Use PETSC_CURRENT to leave any of the values unchanged. Use PETSC_DETERMINE | ||
699 | to set m to the leading dimension and p to the number of columns of B. | ||
700 | |||
701 | Level: intermediate | ||
702 | |||
703 | .seealso: DSGSVDGetDimensions(), DSSetDimensions() | ||
704 | @*/ | ||
705 | 6331 | PetscErrorCode DSGSVDSetDimensions(DS ds,PetscInt m,PetscInt p) | |
706 | { | ||
707 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
6331 | PetscFunctionBegin; |
708 |
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.
|
6331 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
709 |
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.
|
6331 | PetscValidLogicalCollectiveInt(ds,m,2); |
710 |
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.
|
6331 | PetscValidLogicalCollectiveInt(ds,p,3); |
711 |
8/14✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 10 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.
|
6331 | PetscTryMethod(ds,"DSGSVDSetDimensions_C",(DS,PetscInt,PetscInt),(ds,m,p)); |
712 |
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.
|
6331 | PetscFunctionReturn(PETSC_SUCCESS); |
713 | } | ||
714 | |||
715 | 120 | static PetscErrorCode DSGSVDGetDimensions_GSVD(DS ds,PetscInt *m,PetscInt *p) | |
716 | { | ||
717 | 120 | DS_GSVD *ctx = (DS_GSVD*)ds->data; | |
718 | |||
719 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
120 | PetscFunctionBegin; |
720 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
120 | if (m) *m = ctx->m; |
721 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
120 | if (p) *p = ctx->p; |
722 |
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.
|
120 | PetscFunctionReturn(PETSC_SUCCESS); |
723 | } | ||
724 | |||
725 | /*@ | ||
726 | DSGSVDGetDimensions - Returns the number of columns and rows for a DSGSVD. | ||
727 | |||
728 | Not Collective | ||
729 | |||
730 | Input Parameter: | ||
731 | . ds - the direct solver context | ||
732 | |||
733 | Output Parameters: | ||
734 | + m - the number of columns | ||
735 | - p - the number of rows for the second matrix (B) | ||
736 | |||
737 | Level: intermediate | ||
738 | |||
739 | .seealso: DSGSVDSetDimensions() | ||
740 | @*/ | ||
741 | 120 | PetscErrorCode DSGSVDGetDimensions(DS ds,PetscInt *m,PetscInt *p) | |
742 | { | ||
743 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
120 | PetscFunctionBegin; |
744 |
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.
|
120 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
745 |
9/16✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 10 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.
|
120 | PetscUseMethod(ds,"DSGSVDGetDimensions_C",(DS,PetscInt*,PetscInt*),(ds,m,p)); |
746 |
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.
|
120 | PetscFunctionReturn(PETSC_SUCCESS); |
747 | } | ||
748 | |||
749 | 684 | static PetscErrorCode DSDestroy_GSVD(DS ds) | |
750 | { | ||
751 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
684 | PetscFunctionBegin; |
752 |
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.
|
684 | PetscCall(PetscFree(ds->data)); |
753 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
684 | PetscCall(PetscObjectComposeFunction((PetscObject)ds,"DSGSVDSetDimensions_C",NULL)); |
754 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
684 | PetscCall(PetscObjectComposeFunction((PetscObject)ds,"DSGSVDGetDimensions_C",NULL)); |
755 |
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.
|
140 | PetscFunctionReturn(PETSC_SUCCESS); |
756 | } | ||
757 | |||
758 | 10 | static PetscErrorCode DSReallocate_GSVD(DS ds,PetscInt ld) | |
759 | { | ||
760 | 10 | PetscInt i,*perm=ds->perm; | |
761 | |||
762 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
763 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
230 | for (i=0;i<DS_NUM_MAT;i++) { |
764 |
10/12✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
✓ Branch 4 taken 5 times.
✓ Branch 5 taken 5 times.
✓ Branch 6 taken 1 times.
✓ Branch 7 taken 4 times.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
|
220 | if (i!=DS_MAT_A && i!=DS_MAT_B && i!=DS_MAT_X && i!=DS_MAT_U && i!=DS_MAT_V && i!=DS_MAT_T && i!=DS_MAT_D) PetscCall(MatDestroy(&ds->omat[i])); |
765 | } | ||
766 | |||
767 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
10 | PetscCall(DSReallocateMat_Private(ds,DS_MAT_A,ld)); |
768 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
10 | PetscCall(DSReallocateMat_Private(ds,DS_MAT_B,ld)); |
769 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
10 | PetscCall(DSReallocateMat_Private(ds,DS_MAT_X,ld)); |
770 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
10 | PetscCall(DSReallocateMat_Private(ds,DS_MAT_U,ld)); |
771 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
10 | PetscCall(DSReallocateMat_Private(ds,DS_MAT_V,ld)); |
772 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
10 | PetscCall(DSReallocateMat_Private(ds,DS_MAT_T,ld)); |
773 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
10 | PetscCall(DSReallocateMat_Private(ds,DS_MAT_D,ld)); |
774 | |||
775 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
10 | PetscCall(PetscMalloc1(ld,&ds->perm)); |
776 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
10 | PetscCall(PetscArraycpy(ds->perm,perm,ds->ld)); |
777 |
6/8✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
|
10 | PetscCall(PetscFree(perm)); |
778 |
6/12✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
|
2 | PetscFunctionReturn(PETSC_SUCCESS); |
779 | } | ||
780 | |||
781 | /*MC | ||
782 | DSGSVD - Dense Generalized Singular Value Decomposition. | ||
783 | |||
784 | Level: beginner | ||
785 | |||
786 | Notes: | ||
787 | The problem is expressed as A*X = U*C, B*X = V*S, where A and B are | ||
788 | matrices with the same number of columns, m, U and V are orthogonal | ||
789 | (unitary), and X is an mxm invertible matrix. The DS object does not | ||
790 | expose matrices C and S, instead the singular values sigma, which are | ||
791 | the ratios c_i/s_i, are returned in the arguments of DSSolve(). | ||
792 | Note that the number of columns of the returned X, U, V may be smaller | ||
793 | in the case that some c_i or s_i are zero. | ||
794 | |||
795 | The number of rows of A (and U) is the value n passed with DSSetDimensions(). | ||
796 | The number of columns m and the number of rows of B (and V) must be | ||
797 | set via DSGSVDSetDimensions(). | ||
798 | |||
799 | Internally, LAPACK's representation is used, U'*A*Q = C*[0 R], V'*B*Q = S*[0 R], | ||
800 | where X = Q*inv(R) is computed at the end of DSSolve(). | ||
801 | |||
802 | If the compact storage format is selected, then a simplified problem is | ||
803 | solved, where A and B are bidiagonal (possibly with an arrow), and [A;B] | ||
804 | is assumed to have orthonormal columns. We consider two cases: (1) A and B | ||
805 | are square mxm upper bidiagonal, and (2) A is lower bidiagonal with m+1 | ||
806 | rows and B is square upper bidiagonal. In these cases, R=I so it | ||
807 | corresponds to the CS decomposition. The first matrix is stored in two | ||
808 | diagonals of DS_MAT_T, while the second matrix is stored in DS_MAT_D | ||
809 | and the remaining diagonal of DS_MAT_T. | ||
810 | |||
811 | Allowed arguments of DSVectors() are DS_MAT_U, DS_MAT_V and DS_MAT_X. | ||
812 | |||
813 | Used DS matrices: | ||
814 | + DS_MAT_A - first problem matrix | ||
815 | . DS_MAT_B - second problem matrix | ||
816 | . DS_MAT_T - first upper bidiagonal matrix (if compact storage is selected) | ||
817 | . DS_MAT_D - second upper bidiagonal matrix (if compact storage is selected) | ||
818 | . DS_MAT_U - (upper) left generalized singular vectors | ||
819 | . DS_MAT_V - (lower) left generalized singular vectors | ||
820 | - DS_MAT_X - right generalized singular vectors | ||
821 | |||
822 | Implemented methods: | ||
823 | . 0 - Lapack (_ggsvd3 if available, or _ggsvd) | ||
824 | |||
825 | .seealso: DSCreate(), DSSetType(), DSType, DSGSVDSetDimensions() | ||
826 | M*/ | ||
827 | 684 | SLEPC_EXTERN PetscErrorCode DSCreate_GSVD(DS ds) | |
828 | { | ||
829 | 684 | DS_GSVD *ctx; | |
830 | |||
831 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
684 | PetscFunctionBegin; |
832 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
684 | PetscCall(PetscNew(&ctx)); |
833 | 684 | ds->data = (void*)ctx; | |
834 | |||
835 | 684 | ds->ops->allocate = DSAllocate_GSVD; | |
836 | 684 | ds->ops->view = DSView_GSVD; | |
837 | 684 | ds->ops->vectors = DSVectors_GSVD; | |
838 | 684 | ds->ops->sort = DSSort_GSVD; | |
839 | 684 | ds->ops->solve[0] = DSSolve_GSVD; | |
840 | #if !defined(PETSC_HAVE_MPIUNI) | ||
841 | 684 | ds->ops->synchronize = DSSynchronize_GSVD; | |
842 | #endif | ||
843 | 684 | ds->ops->truncate = DSTruncate_GSVD; | |
844 | 684 | ds->ops->update = DSUpdateExtraRow_GSVD; | |
845 | 684 | ds->ops->cond = DSCond_GSVD; | |
846 | 684 | ds->ops->matgetsize = DSMatGetSize_GSVD; | |
847 | 684 | ds->ops->destroy = DSDestroy_GSVD; | |
848 | 684 | ds->ops->reallocate = DSReallocate_GSVD; | |
849 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
684 | PetscCall(PetscObjectComposeFunction((PetscObject)ds,"DSGSVDSetDimensions_C",DSGSVDSetDimensions_GSVD)); |
850 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
684 | PetscCall(PetscObjectComposeFunction((PetscObject)ds,"DSGSVDGetDimensions_C",DSGSVDGetDimensions_GSVD)); |
851 |
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.
|
140 | PetscFunctionReturn(PETSC_SUCCESS); |
852 | } | ||
853 |