Line | Branch | Exec | Source |
---|---|---|---|
1 | /* | ||
2 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
3 | SLEPc - Scalable Library for Eigenvalue Problem Computations | ||
4 | Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain | ||
5 | |||
6 | This file is part of SLEPc. | ||
7 | SLEPc is distributed under a 2-clause BSD license (see LICENSE). | ||
8 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
9 | */ | ||
10 | /* | ||
11 | Private DS routines | ||
12 | */ | ||
13 | |||
14 | #include <slepc/private/dsimpl.h> /*I "slepcds.h" I*/ | ||
15 | #include <slepcblaslapack.h> | ||
16 | |||
17 | 106445 | PetscErrorCode DSAllocateMat_Private(DS ds,DSMatType m) | |
18 | { | ||
19 | 106445 | PetscInt n,d,rows=0,cols; | |
20 | 106445 | PetscBool ispep,isnep; | |
21 | |||
22 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
106445 | PetscFunctionBegin; |
23 | 106445 | n = ds->ld; | |
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.
|
106445 | PetscCall(PetscObjectTypeCompare((PetscObject)ds,DSPEP,&ispep)); |
25 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
106445 | if (ispep) { |
26 |
6/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 2 times.
|
13541 | if (m==DS_MAT_A || m==DS_MAT_B || m==DS_MAT_W || m==DS_MAT_U || m==DS_MAT_X || m==DS_MAT_Y) { |
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.
|
12906 | PetscCall(DSPEPGetDegree(ds,&d)); |
28 | 12906 | n = d*ds->ld; | |
29 | } | ||
30 | } else { | ||
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.
|
92904 | PetscCall(PetscObjectTypeCompare((PetscObject)ds,DSNEP,&isnep)); |
32 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
92904 | if (isnep) { |
33 |
5/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 10 times.
|
7033 | if (m==DS_MAT_Q || m==DS_MAT_Z || m==DS_MAT_U || m==DS_MAT_V || m==DS_MAT_X || m==DS_MAT_Y) { |
34 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
849 | PetscCall(DSNEPGetMinimality(ds,&d)); |
35 | 849 | n = d*ds->ld; | |
36 | } | ||
37 | } | ||
38 | } | ||
39 | 106445 | cols = n; | |
40 | |||
41 |
5/5✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 10 times.
|
106445 | switch (m) { |
42 | case DS_MAT_T: | ||
43 | cols = PetscDefined(USE_COMPLEX)? 2: 3; | ||
44 | rows = n; | ||
45 | break; | ||
46 | 1437 | case DS_MAT_D: | |
47 | 1437 | cols = 1; | |
48 | 1437 | rows = n; | |
49 | 1437 | break; | |
50 | 19444 | case DS_MAT_X: | |
51 | 19444 | rows = ds->ld; | |
52 | 19444 | break; | |
53 | 912 | case DS_MAT_Y: | |
54 | 912 | rows = ds->ld; | |
55 | 912 | break; | |
56 | 78656 | default: | |
57 | 78656 | rows = n; | |
58 | } | ||
59 |
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.
|
106445 | if (ds->omat[m]) PetscCall(MatZeroEntries(ds->omat[m])); |
60 | else { | ||
61 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
45813 | PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,rows,cols,NULL,&ds->omat[m])); |
62 | } | ||
63 |
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.
|
21146 | PetscFunctionReturn(PETSC_SUCCESS); |
64 | } | ||
65 | |||
66 | 331 | PetscErrorCode DSReallocateMat_Private(DS ds,DSMatType mat,PetscInt ld) | |
67 | { | ||
68 | 331 | Mat M; | |
69 | 331 | PetscInt i,m,n,rows,cols; | |
70 | 331 | PetscScalar *dst; | |
71 | 331 | PetscReal *rdst; | |
72 | 331 | const PetscScalar *src; | |
73 | 331 | const PetscReal *rsrc; | |
74 | |||
75 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
331 | PetscFunctionBegin; |
76 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
331 | PetscCall(MatGetSize(ds->omat[mat],&m,&n)); |
77 | 331 | rows = ld; | |
78 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
331 | cols = m==n? ld: n; |
79 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
331 | PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,rows,cols,NULL,&M)); |
80 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
331 | PetscCall(MatDenseGetArrayRead(ds->omat[mat],&src)); |
81 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
331 | PetscCall(MatDenseGetArray(M,&dst)); |
82 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
126 | if (mat==DS_MAT_T && PetscDefined(USE_COMPLEX)) { |
83 | 47 | rsrc = (const PetscReal*)src; | |
84 | 47 | rdst = (PetscReal*)dst; | |
85 |
7/8✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 4 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
✓ Branch 6 taken 1 times.
✓ Branch 7 taken 1 times.
|
188 | for (i=0;i<3;i++) PetscCall(PetscArraycpy(rdst+i*ld,rsrc+i*ds->ld,ds->ld)); |
86 | } else { | ||
87 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
2973 | for (i=0;i<n;i++) PetscCall(PetscArraycpy(dst+i*ld,src+i*ds->ld,ds->ld)); |
88 | } | ||
89 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
331 | PetscCall(MatDenseRestoreArrayRead(ds->omat[mat],&src)); |
90 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
331 | PetscCall(MatDenseRestoreArray(M,&dst)); |
91 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
331 | PetscCall(MatDestroy(&ds->omat[mat])); |
92 | 331 | ds->omat[mat] = M; | |
93 |
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.
|
331 | PetscFunctionReturn(PETSC_SUCCESS); |
94 | } | ||
95 | |||
96 | 834039 | PetscErrorCode DSAllocateWork_Private(DS ds,PetscInt s,PetscInt r,PetscInt i) | |
97 | { | ||
98 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
834039 | PetscFunctionBegin; |
99 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
834039 | if (s>ds->lwork) { |
100 |
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.
|
16375 | PetscCall(PetscFree(ds->work)); |
101 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
16375 | PetscCall(PetscMalloc1(s,&ds->work)); |
102 | 16375 | ds->lwork = s; | |
103 | } | ||
104 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
834039 | if (r>ds->lrwork) { |
105 |
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.
|
10180 | PetscCall(PetscFree(ds->rwork)); |
106 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10180 | PetscCall(PetscMalloc1(r,&ds->rwork)); |
107 | 10180 | ds->lrwork = r; | |
108 | } | ||
109 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
834039 | if (i>ds->liwork) { |
110 |
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.
|
6629 | PetscCall(PetscFree(ds->iwork)); |
111 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6629 | PetscCall(PetscMalloc1(i,&ds->iwork)); |
112 | 6629 | ds->liwork = i; | |
113 | } | ||
114 |
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.
|
164537 | PetscFunctionReturn(PETSC_SUCCESS); |
115 | } | ||
116 | |||
117 | /*@ | ||
118 | DSViewMat - Prints one of the internal DS matrices. | ||
119 | |||
120 | Collective | ||
121 | |||
122 | Input Parameters: | ||
123 | + ds - the direct solver context | ||
124 | . viewer - visualization context | ||
125 | - m - matrix to display | ||
126 | |||
127 | Note: | ||
128 | Works only for ascii viewers. Set the viewer in Matlab format if | ||
129 | want to paste into Matlab. | ||
130 | |||
131 | Level: developer | ||
132 | |||
133 | .seealso: DSView() | ||
134 | @*/ | ||
135 | 300 | PetscErrorCode DSViewMat(DS ds,PetscViewer viewer,DSMatType m) | |
136 | { | ||
137 | 300 | PetscInt i,j,rows,cols; | |
138 | 300 | const PetscScalar *M=NULL,*v; | |
139 | 300 | PetscViewerFormat format; | |
140 | #if defined(PETSC_USE_COMPLEX) | ||
141 | 150 | PetscBool allreal = PETSC_TRUE; | |
142 | 150 | const PetscReal *vr; | |
143 | #endif | ||
144 | |||
145 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
300 | PetscFunctionBegin; |
146 |
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.
|
300 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
147 |
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.
|
300 | PetscValidLogicalCollectiveEnum(ds,m,3); |
148 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
300 | DSCheckValidMat(ds,m,3); |
149 |
1/8✗ 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.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
300 | if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ds),&viewer)); |
150 |
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.
|
300 | PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); |
151 |
13/32✓ 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 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.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
|
300 | PetscCheckSameComm(ds,1,viewer,2); |
152 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
300 | PetscCall(PetscViewerGetFormat(viewer,&format)); |
153 |
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 || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(PETSC_SUCCESS); |
154 |
4/6✓ Branch 0 taken 2 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(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); |
155 |
4/6✓ Branch 0 taken 2 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(DSMatGetSize(ds,m,&rows,&cols)); |
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.
|
300 | PetscCall(MatDenseGetArrayRead(ds->omat[m],&M)); |
157 | #if defined(PETSC_USE_COMPLEX) | ||
158 | /* determine if matrix has all real values */ | ||
159 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
150 | if (m!=DS_MAT_T && m!=DS_MAT_D) { |
160 | /* determine if matrix has all real values */ | ||
161 | 145 | v = M; | |
162 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
2020 | for (i=0;i<rows;i++) |
163 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
20775 | for (j=0;j<cols;j++) |
164 |
1/2✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
|
18900 | if (PetscImaginaryPart(v[i+j*ds->ld])) { allreal = PETSC_FALSE; break; } |
165 | } | ||
166 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
150 | if (m==DS_MAT_T) cols=3; |
167 | #endif | ||
168 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
300 | if (format == PETSC_VIEWER_ASCII_MATLAB) { |
169 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer,"%% Size = %" PetscInt_FMT " %" PetscInt_FMT "\n",rows,cols)); | |
170 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer,"%s = [\n",DSMatName[m])); | |
171 |
4/6✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✓ Branch 5 taken 1 times.
|
300 | } else PetscCall(PetscViewerASCIIPrintf(viewer,"Matrix %s =\n",DSMatName[m])); |
172 | |||
173 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
4200 | for (i=0;i<rows;i++) { |
174 | 3900 | v = M+i; | |
175 | #if defined(PETSC_USE_COMPLEX) | ||
176 | 1950 | vr = (const PetscReal*)M+i; /* handle compact storage, 2nd column is in imaginary part of 1st column */ | |
177 | #endif | ||
178 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
41850 | for (j=0;j<cols;j++) { |
179 | #if defined(PETSC_USE_COMPLEX) | ||
180 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
18975 | if (m!=DS_MAT_T && m!=DS_MAT_D) { |
181 |
5/8✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
|
18900 | if (allreal) PetscCall(PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(*v))); |
182 | ✗ | else PetscCall(PetscViewerASCIIPrintf(viewer,"%18.16e%+18.16ei ",(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v))); | |
183 |
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.
|
75 | } else PetscCall(PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)*vr)); |
184 | #else | ||
185 |
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.
|
18975 | PetscCall(PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)*v)); |
186 | #endif | ||
187 | 37950 | v += ds->ld; | |
188 | #if defined(PETSC_USE_COMPLEX) | ||
189 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
18975 | if (m==DS_MAT_T) vr += ds->ld; |
190 | #endif | ||
191 | } | ||
192 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3900 | PetscCall(PetscViewerASCIIPrintf(viewer,"\n")); |
193 | } | ||
194 |
4/6✓ Branch 0 taken 2 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(MatDenseRestoreArrayRead(ds->omat[m],&M)); |
195 | |||
196 |
1/8✗ 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.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
300 | if (format == PETSC_VIEWER_ASCII_MATLAB) PetscCall(PetscViewerASCIIPrintf(viewer,"];\n")); |
197 |
4/6✓ Branch 0 taken 2 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(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); |
198 |
4/6✓ Branch 0 taken 2 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(PetscViewerFlush(viewer)); |
199 |
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); |
200 | } | ||
201 | |||
202 | 191158 | PetscErrorCode DSSortEigenvalues_Private(DS ds,PetscScalar *wr,PetscScalar *wi,PetscInt *perm,PetscBool isghiep) | |
203 | { | ||
204 | 191158 | PetscScalar re,im,wi0; | |
205 | 191158 | PetscInt n,i,j,result,tmp1,tmp2=0,d=1; | |
206 | |||
207 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
191158 | PetscFunctionBegin; |
208 | 191158 | n = ds->t; /* sort only first t pairs if truncated */ | |
209 | /* insertion sort */ | ||
210 | 191158 | i=ds->l+1; | |
211 | #if !defined(PETSC_USE_COMPLEX) | ||
212 |
4/4✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
|
94009 | if (wi && wi[perm[i-1]]!=0.0) i++; /* initial value is complex */ |
213 | #else | ||
214 |
4/4✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
|
97149 | if (isghiep && PetscImaginaryPart(wr[perm[i-1]])!=0.0) i++; |
215 | #endif | ||
216 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2355411 | for (;i<n;i+=d) { |
217 | 2164253 | re = wr[perm[i]]; | |
218 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2164253 | if (wi) im = wi[perm[i]]; |
219 | else im = 0.0; | ||
220 | 2164253 | tmp1 = perm[i]; | |
221 | #if !defined(PETSC_USE_COMPLEX) | ||
222 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
1126429 | if (im!=0.0) { d = 2; tmp2 = perm[i+1]; } |
223 | else d = 1; | ||
224 | #else | ||
225 |
4/4✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
|
1037824 | if (isghiep && PetscImaginaryPart(re)!=0.0) { d = 2; tmp2 = perm[i+1]; } |
226 | else d = 1; | ||
227 | #endif | ||
228 | 2164253 | j = i-1; | |
229 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2164253 | if (wi) wi0 = wi[perm[j]]; |
230 | else wi0 = 0.0; | ||
231 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2164253 | PetscCall(SlepcSCCompare(ds->sc,re,im,wr[perm[j]],wi0,&result)); |
232 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
8331444 | while (result<0 && j>=ds->l) { |
233 | 7710433 | perm[j+d] = perm[j]; | |
234 | 7710433 | j--; | |
235 | #if !defined(PETSC_USE_COMPLEX) | ||
236 |
4/4✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
|
4895384 | if (wi && wi[perm[j+1]]!=0) |
237 | #else | ||
238 |
4/4✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
|
2815049 | if (isghiep && PetscImaginaryPart(wr[perm[j+1]])!=0) |
239 | #endif | ||
240 | 205427 | { perm[j+d] = perm[j]; j--; } | |
241 | |||
242 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
7710433 | if (j>=ds->l) { |
243 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
7089422 | if (wi) wi0 = wi[perm[j]]; |
244 | else wi0 = 0.0; | ||
245 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
16964108 | PetscCall(SlepcSCCompare(ds->sc,re,im,wr[perm[j]],wi0,&result)); |
246 | } | ||
247 | } | ||
248 | 2164253 | perm[j+1] = tmp1; | |
249 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2164253 | if (d==2) perm[j+2] = tmp2; |
250 | } | ||
251 |
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.
|
37403 | PetscFunctionReturn(PETSC_SUCCESS); |
252 | } | ||
253 | |||
254 | 68473 | PetscErrorCode DSSortEigenvaluesReal_Private(DS ds,PetscReal *eig,PetscInt *perm) | |
255 | { | ||
256 | 68473 | PetscScalar re; | |
257 | 68473 | PetscInt i,j,result,tmp,l,n; | |
258 | |||
259 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
68473 | PetscFunctionBegin; |
260 | 68473 | n = ds->t; /* sort only first t pairs if truncated */ | |
261 | 68473 | l = ds->l; | |
262 | /* insertion sort */ | ||
263 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
942955 | for (i=l+1;i<n;i++) { |
264 | 874482 | re = eig[perm[i]]; | |
265 | 874482 | j = i-1; | |
266 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
874482 | PetscCall(SlepcSCCompare(ds->sc,re,0.0,eig[perm[j]],0.0,&result)); |
267 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
5266828 | while (result<0 && j>=l) { |
268 | 4934054 | tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--; | |
269 |
9/10✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 8 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✓ Branch 9 taken 2 times.
|
10742590 | if (j>=l) PetscCall(SlepcSCCompare(ds->sc,re,0.0,eig[perm[j]],0.0,&result)); |
270 | } | ||
271 | } | ||
272 |
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.
|
12639 | PetscFunctionReturn(PETSC_SUCCESS); |
273 | } | ||
274 | |||
275 | /* | ||
276 | Permute comumns [istart..iend-1] of [mat] according to perm. Columns have length n | ||
277 | */ | ||
278 | 218340 | PetscErrorCode DSPermuteColumns_Private(DS ds,PetscInt istart,PetscInt iend,PetscInt n,DSMatType mat,PetscInt *perm) | |
279 | { | ||
280 | 218340 | PetscInt i,j,k,p,ld; | |
281 | 218340 | PetscScalar *Q,rtmp; | |
282 | |||
283 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
218340 | PetscFunctionBegin; |
284 | 218340 | ld = ds->ld; | |
285 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
218340 | PetscCall(MatDenseGetArray(ds->omat[mat],&Q)); |
286 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2995746 | for (i=istart;i<iend;i++) { |
287 | 2777406 | p = perm[i]; | |
288 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2777406 | if (p != i) { |
289 | 980633 | j = i + 1; | |
290 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
8115411 | while (perm[j] != i) j++; |
291 | 980633 | perm[j] = p; perm[i] = i; | |
292 | /* swap columns i and j */ | ||
293 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
26998942 | for (k=0;k<n;k++) { |
294 | 26018309 | rtmp = Q[k+p*ld]; Q[k+p*ld] = Q[k+i*ld]; Q[k+i*ld] = rtmp; | |
295 | } | ||
296 | } | ||
297 | } | ||
298 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
218340 | PetscCall(MatDenseRestoreArray(ds->omat[mat],&Q)); |
299 |
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.
|
42447 | PetscFunctionReturn(PETSC_SUCCESS); |
300 | } | ||
301 | |||
302 | /* | ||
303 | The same as DSPermuteColumns_Private but for two matrices [mat1] and [mat2] | ||
304 | */ | ||
305 | 3134 | PetscErrorCode DSPermuteColumnsTwo_Private(DS ds,PetscInt istart,PetscInt iend,PetscInt n,DSMatType mat1,DSMatType mat2,PetscInt *perm) | |
306 | { | ||
307 | 3134 | PetscInt i,j,k,p,ld; | |
308 | 3134 | PetscScalar *Q,*Z,rtmp,rtmp2; | |
309 | |||
310 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
3134 | PetscFunctionBegin; |
311 | 3134 | ld = ds->ld; | |
312 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3134 | PetscCall(MatDenseGetArray(ds->omat[mat1],&Q)); |
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.
|
3134 | PetscCall(MatDenseGetArray(ds->omat[mat2],&Z)); |
314 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
76187 | for (i=istart;i<iend;i++) { |
315 | 73053 | p = perm[i]; | |
316 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
73053 | if (p != i) { |
317 | 52245 | j = i + 1; | |
318 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
516273 | while (perm[j] != i) j++; |
319 | 52245 | perm[j] = p; perm[i] = i; | |
320 | /* swap columns i and j */ | ||
321 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
853739 | for (k=0;k<n;k++) { |
322 | 801494 | rtmp = Q[k+p*ld]; Q[k+p*ld] = Q[k+i*ld]; Q[k+i*ld] = rtmp; | |
323 | 801494 | rtmp2 = Z[k+p*ld]; Z[k+p*ld] = Z[k+i*ld]; Z[k+i*ld] = rtmp2; | |
324 | } | ||
325 | } | ||
326 | } | ||
327 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3134 | PetscCall(MatDenseRestoreArray(ds->omat[mat1],&Q)); |
328 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3134 | PetscCall(MatDenseRestoreArray(ds->omat[mat2],&Z)); |
329 |
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.
|
659 | PetscFunctionReturn(PETSC_SUCCESS); |
330 | } | ||
331 | |||
332 | /* | ||
333 | Permute rows [istart..iend-1] of [mat] according to perm. Rows have length m | ||
334 | */ | ||
335 | ✗ | PetscErrorCode DSPermuteRows_Private(DS ds,PetscInt istart,PetscInt iend,PetscInt m,DSMatType mat,PetscInt *perm) | |
336 | { | ||
337 | ✗ | PetscInt i,j,k,p,ld; | |
338 | ✗ | PetscScalar *Q,rtmp; | |
339 | |||
340 | ✗ | PetscFunctionBegin; | |
341 | ✗ | ld = ds->ld; | |
342 | ✗ | PetscCall(MatDenseGetArray(ds->omat[mat],&Q)); | |
343 | ✗ | for (i=istart;i<iend;i++) { | |
344 | ✗ | p = perm[i]; | |
345 | ✗ | if (p != i) { | |
346 | ✗ | j = i + 1; | |
347 | ✗ | while (perm[j] != i) j++; | |
348 | ✗ | perm[j] = p; perm[i] = i; | |
349 | /* swap rows i and j */ | ||
350 | ✗ | for (k=0;k<m;k++) { | |
351 | ✗ | rtmp = Q[p+k*ld]; Q[p+k*ld] = Q[i+k*ld]; Q[i+k*ld] = rtmp; | |
352 | } | ||
353 | } | ||
354 | } | ||
355 | ✗ | PetscCall(MatDenseRestoreArray(ds->omat[mat],&Q)); | |
356 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
357 | } | ||
358 | |||
359 | /* | ||
360 | Permute columns [istart..iend-1] of [mat1] and [mat2] according to perm. | ||
361 | Columns of [mat1] have length n, columns of [mat2] have length m | ||
362 | */ | ||
363 | 25134 | PetscErrorCode DSPermuteBoth_Private(DS ds,PetscInt istart,PetscInt iend,PetscInt n,PetscInt m,DSMatType mat1,DSMatType mat2,PetscInt *perm) | |
364 | { | ||
365 | 25134 | PetscInt i,j,k,p,ld; | |
366 | 25134 | PetscScalar *U,*V,rtmp; | |
367 | |||
368 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
25134 | PetscFunctionBegin; |
369 | 25134 | ld = ds->ld; | |
370 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
25134 | PetscCall(MatDenseGetArray(ds->omat[mat1],&U)); |
371 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
25134 | PetscCall(MatDenseGetArray(ds->omat[mat2],&V)); |
372 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
304047 | for (i=istart;i<iend;i++) { |
373 | 278913 | p = perm[i]; | |
374 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
278913 | if (p != i) { |
375 | 62754 | j = i + 1; | |
376 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
618372 | while (perm[j] != i) j++; |
377 | 62754 | perm[j] = p; perm[i] = i; | |
378 | /* swap columns i and j of U */ | ||
379 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
3271322 | for (k=0;k<n;k++) { |
380 | 3208568 | rtmp = U[k+p*ld]; U[k+p*ld] = U[k+i*ld]; U[k+i*ld] = rtmp; | |
381 | } | ||
382 | /* swap columns i and j of V */ | ||
383 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2058109 | for (k=0;k<m;k++) { |
384 | 1995355 | rtmp = V[k+p*ld]; V[k+p*ld] = V[k+i*ld]; V[k+i*ld] = rtmp; | |
385 | } | ||
386 | } | ||
387 | } | ||
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.
|
25134 | PetscCall(MatDenseRestoreArray(ds->omat[mat1],&U)); |
389 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
25134 | PetscCall(MatDenseRestoreArray(ds->omat[mat2],&V)); |
390 |
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.
|
4187 | PetscFunctionReturn(PETSC_SUCCESS); |
391 | } | ||
392 | |||
393 | /*@ | ||
394 | DSSetIdentity - Copy the identity (a diagonal matrix with ones) on the | ||
395 | active part of a matrix. | ||
396 | |||
397 | Logically Collective | ||
398 | |||
399 | Input Parameters: | ||
400 | + ds - the direct solver context | ||
401 | - mat - the matrix to modify | ||
402 | |||
403 | Level: intermediate | ||
404 | |||
405 | .seealso: DSGetMat() | ||
406 | @*/ | ||
407 | 106275 | PetscErrorCode DSSetIdentity(DS ds,DSMatType mat) | |
408 | { | ||
409 | 106275 | PetscScalar *A; | |
410 | 106275 | PetscInt i,ld,n,l; | |
411 | |||
412 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
106275 | PetscFunctionBegin; |
413 |
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.
|
106275 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
414 |
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.
|
106275 | PetscValidLogicalCollectiveEnum(ds,mat,2); |
415 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
106275 | DSCheckValidMat(ds,mat,2); |
416 | |||
417 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
106275 | PetscCall(DSGetDimensions(ds,&n,&l,NULL,NULL)); |
418 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
106275 | PetscCall(DSGetLeadingDimension(ds,&ld)); |
419 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
106275 | PetscCall(PetscLogEventBegin(DS_Other,ds,0,0,0)); |
420 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
106275 | PetscCall(MatDenseGetArray(ds->omat[mat],&A)); |
421 |
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.
|
106275 | PetscCall(PetscArrayzero(A+l*ld,ld*(n-l))); |
422 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1479905 | for (i=l;i<n;i++) A[i+i*ld] = 1.0; |
423 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
106275 | PetscCall(MatDenseRestoreArray(ds->omat[mat],&A)); |
424 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
106275 | PetscCall(PetscLogEventEnd(DS_Other,ds,0,0,0)); |
425 |
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.
|
20517 | PetscFunctionReturn(PETSC_SUCCESS); |
426 | } | ||
427 | |||
428 | /*@ | ||
429 | DSOrthogonalize - Orthogonalize the columns of a matrix. | ||
430 | |||
431 | Logically Collective | ||
432 | |||
433 | Input Parameters: | ||
434 | + ds - the direct solver context | ||
435 | . mat - a matrix | ||
436 | - cols - number of columns to orthogonalize (starting from column zero) | ||
437 | |||
438 | Output Parameter: | ||
439 | . lindcols - (optional) number of linearly independent columns of the matrix | ||
440 | |||
441 | Level: developer | ||
442 | |||
443 | .seealso: DSPseudoOrthogonalize() | ||
444 | @*/ | ||
445 | 16907 | PetscErrorCode DSOrthogonalize(DS ds,DSMatType mat,PetscInt cols,PetscInt *lindcols) | |
446 | { | ||
447 | 16907 | PetscInt n,l,ld; | |
448 | 16907 | PetscBLASInt ld_,rA,cA,info,ltau,lw; | |
449 | 16907 | PetscScalar *A,*tau,*w,saux,dummy; | |
450 | |||
451 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
16907 | PetscFunctionBegin; |
452 |
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.
|
16907 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
453 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
16907 | DSCheckAlloc(ds,1); |
454 |
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.
|
16907 | PetscValidLogicalCollectiveEnum(ds,mat,2); |
455 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
16907 | DSCheckValidMat(ds,mat,2); |
456 |
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.
|
16907 | PetscValidLogicalCollectiveInt(ds,cols,3); |
457 | |||
458 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
16907 | PetscCall(DSGetDimensions(ds,&n,&l,NULL,NULL)); |
459 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
16907 | PetscCall(DSGetLeadingDimension(ds,&ld)); |
460 | 16907 | n = n - l; | |
461 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
16907 | PetscCheck(cols<=n,PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid number of columns"); |
462 |
3/16✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 2 times.
✗ 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.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
16907 | if (n == 0 || cols == 0) PetscFunctionReturn(PETSC_SUCCESS); |
463 | |||
464 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
16907 | PetscCall(PetscLogEventBegin(DS_Other,ds,0,0,0)); |
465 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
16907 | PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); |
466 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
16907 | PetscCall(MatDenseGetArray(ds->omat[mat],&A)); |
467 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
16907 | PetscCall(PetscBLASIntCast(PetscMin(cols,n),<au)); |
468 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
16907 | PetscCall(PetscBLASIntCast(ld,&ld_)); |
469 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
16907 | PetscCall(PetscBLASIntCast(n,&rA)); |
470 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
16907 | PetscCall(PetscBLASIntCast(cols,&cA)); |
471 | 16907 | lw = -1; | |
472 |
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.
|
16907 | PetscCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&rA,&cA,A,&ld_,&dummy,&saux,&lw,&info)); |
473 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
16907 | SlepcCheckLapackInfo("geqrf",info); |
474 | 16907 | lw = (PetscBLASInt)PetscRealPart(saux); | |
475 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
16907 | PetscCall(DSAllocateWork_Private(ds,lw+ltau,0,0)); |
476 | 16907 | tau = ds->work; | |
477 | 16907 | w = &tau[ltau]; | |
478 |
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.
|
16907 | PetscCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&rA,&cA,&A[ld*l+l],&ld_,tau,w,&lw,&info)); |
479 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
16907 | SlepcCheckLapackInfo("geqrf",info); |
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.
|
16907 | PetscCallBLAS("LAPACKorgqr",LAPACKorgqr_(&rA,<au,<au,&A[ld*l+l],&ld_,tau,w,&lw,&info)); |
481 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
16907 | SlepcCheckLapackInfo("orgqr",info); |
482 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
16907 | if (lindcols) *lindcols = ltau; |
483 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
16907 | PetscCall(PetscFPTrapPop()); |
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.
|
16907 | PetscCall(MatDenseRestoreArray(ds->omat[mat],&A)); |
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.
|
16907 | PetscCall(PetscLogEventEnd(DS_Other,ds,0,0,0)); |
486 |
2/4✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
|
16907 | PetscCall(PetscObjectStateIncrease((PetscObject)ds)); |
487 |
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.
|
16907 | PetscFunctionReturn(PETSC_SUCCESS); |
488 | } | ||
489 | |||
490 | /* | ||
491 | Compute C <- a*A*B + b*C, where | ||
492 | ldC, the leading dimension of C, | ||
493 | ldA, the leading dimension of A, | ||
494 | rA, cA, rows and columns of A, | ||
495 | At, if true use the transpose of A instead, | ||
496 | ldB, the leading dimension of B, | ||
497 | rB, cB, rows and columns of B, | ||
498 | Bt, if true use the transpose of B instead | ||
499 | */ | ||
500 | 550 | static PetscErrorCode SlepcMatDenseMult(PetscScalar *C,PetscInt _ldC,PetscScalar b,PetscScalar a,const PetscScalar *A,PetscInt _ldA,PetscInt rA,PetscInt cA,PetscBool At,const PetscScalar *B,PetscInt _ldB,PetscInt rB,PetscInt cB,PetscBool Bt) | |
501 | { | ||
502 | 550 | PetscInt tmp; | |
503 | 550 | PetscBLASInt m, n, k, ldA, ldB, ldC; | |
504 | 550 | const char *N = "N", *T = "C", *qA = N, *qB = N; | |
505 | |||
506 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
550 | PetscFunctionBegin; |
507 |
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.
|
550 | if ((rA == 0) || (cB == 0)) PetscFunctionReturn(PETSC_SUCCESS); |
508 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
550 | PetscAssertPointer(C,1); |
509 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
550 | PetscAssertPointer(A,5); |
510 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
550 | PetscAssertPointer(B,10); |
511 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
550 | PetscCall(PetscBLASIntCast(_ldA,&ldA)); |
512 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
550 | PetscCall(PetscBLASIntCast(_ldB,&ldB)); |
513 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
550 | PetscCall(PetscBLASIntCast(_ldC,&ldC)); |
514 | |||
515 | /* Transpose if needed */ | ||
516 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
550 | if (At) tmp = rA, rA = cA, cA = tmp, qA = T; |
517 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
550 | if (Bt) tmp = rB, rB = cB, cB = tmp, qB = T; |
518 | |||
519 | /* Check size */ | ||
520 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
550 | PetscCheck(cA==rB,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Matrix dimensions do not match"); |
521 | |||
522 | /* Do stub */ | ||
523 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
550 | if ((rA == 1) && (cA == 1) && (cB == 1)) { |
524 | ✗ | if (!At && !Bt) *C = *A * *B; | |
525 | ✗ | else if (At && !Bt) *C = PetscConj(*A) * *B; | |
526 | ✗ | else if (!At && Bt) *C = *A * PetscConj(*B); | |
527 | ✗ | else *C = PetscConj(*A) * PetscConj(*B); | |
528 | ✗ | m = n = k = 1; | |
529 | } else { | ||
530 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
550 | PetscCall(PetscBLASIntCast(rA,&m)); |
531 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
550 | PetscCall(PetscBLASIntCast(cB,&n)); |
532 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
550 | PetscCall(PetscBLASIntCast(cA,&k)); |
533 |
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.
|
550 | PetscCallBLAS("BLASgemm",BLASgemm_(qA,qB,&m,&n,&k,&a,(PetscScalar*)A,&ldA,(PetscScalar*)B,&ldB,&b,C,&ldC)); |
534 | } | ||
535 | |||
536 |
3/6✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
550 | PetscCall(PetscLogFlops(2.0*m*n*k)); |
537 |
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.
|
550 | PetscFunctionReturn(PETSC_SUCCESS); |
538 | } | ||
539 | |||
540 | /*@ | ||
541 | DSPseudoOrthogonalize - Orthogonalize the columns of a matrix with Modified | ||
542 | Gram-Schmidt in an indefinite inner product space defined by a signature. | ||
543 | |||
544 | Logically Collective | ||
545 | |||
546 | Input Parameters: | ||
547 | + ds - the direct solver context | ||
548 | . mat - the matrix | ||
549 | . cols - number of columns to orthogonalize (starting from column zero) | ||
550 | - s - the signature that defines the inner product | ||
551 | |||
552 | Output Parameters: | ||
553 | + lindcols - (optional) linearly independent columns of the matrix | ||
554 | - ns - (optional) the new signature of the vectors | ||
555 | |||
556 | Note: | ||
557 | After the call the matrix satisfies A'*s*A = ns. | ||
558 | |||
559 | Level: developer | ||
560 | |||
561 | .seealso: DSOrthogonalize() | ||
562 | @*/ | ||
563 | 10 | PetscErrorCode DSPseudoOrthogonalize(DS ds,DSMatType mat,PetscInt cols,PetscReal s[],PetscInt *lindcols,PetscReal ns[]) | |
564 | { | ||
565 | 10 | PetscInt i,j,k,l,n,ld; | |
566 | 10 | PetscBLASInt info,one=1,zero=0,rA_,ld_; | |
567 | 10 | PetscScalar *A,*A_,*m,*h,nr0; | |
568 | 10 | PetscReal nr_o,nr,nr_abs,*ns_,done=1.0; | |
569 | |||
570 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
571 |
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.
|
10 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
572 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
10 | DSCheckAlloc(ds,1); |
573 |
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.
|
10 | PetscValidLogicalCollectiveEnum(ds,mat,2); |
574 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
10 | DSCheckValidMat(ds,mat,2); |
575 |
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.
|
10 | PetscValidLogicalCollectiveInt(ds,cols,3); |
576 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
10 | PetscAssertPointer(s,4); |
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.
|
10 | PetscCall(DSGetDimensions(ds,&n,&l,NULL,NULL)); |
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.
|
10 | PetscCall(DSGetLeadingDimension(ds,&ld)); |
579 | 10 | n = n - l; | |
580 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
10 | PetscCheck(cols<=n,PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid number of columns"); |
581 |
3/16✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 2 times.
✗ 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.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
10 | if (n == 0 || cols == 0) PetscFunctionReturn(PETSC_SUCCESS); |
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.
|
10 | PetscCall(PetscLogEventBegin(DS_Other,ds,0,0,0)); |
583 |
4/6✓ Branch 0 taken 2 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(PetscBLASIntCast(n,&rA_)); |
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.
|
10 | PetscCall(PetscBLASIntCast(ld,&ld_)); |
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.
|
10 | PetscCall(MatDenseGetArray(ds->omat[mat],&A_)); |
586 | 10 | A = &A_[ld*l+l]; | |
587 |
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.
|
20 | PetscCall(DSAllocateWork_Private(ds,n+cols,ns?0:cols,0)); |
588 | 10 | m = ds->work; | |
589 | 10 | h = &m[n]; | |
590 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
10 | ns_ = ns ? ns : ds->rwork; |
591 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
110 | for (i=0; i<cols; i++) { |
592 | /* m <- diag(s)*A[i] */ | ||
593 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1100 | for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld]; |
594 | /* nr_o <- mynorm(A[i]'*m), mynorm(x) = sign(x)*sqrt(|x|) */ | ||
595 |
4/6✓ Branch 0 taken 2 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(SlepcMatDenseMult(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE)); |
596 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
|
100 | nr = nr_o = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0)); |
597 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
160 | for (j=0; j<3 && i>0; j++) { |
598 | /* h <- A[0:i-1]'*m */ | ||
599 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
150 | PetscCall(SlepcMatDenseMult(h,i,0.0,1.0,A,ld,n,i,PETSC_TRUE,m,n,n,1,PETSC_FALSE)); |
600 | /* h <- diag(ns)*h */ | ||
601 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
990 | for (k=0; k<i; k++) h[k] *= ns_[k]; |
602 | /* A[i] <- A[i] - A[0:i-1]*h */ | ||
603 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
150 | PetscCall(SlepcMatDenseMult(&A[ld*i],ld,1.0,-1.0,A,ld,n,i,PETSC_FALSE,h,i,i,1,PETSC_FALSE)); |
604 | /* m <- diag(s)*A[i] */ | ||
605 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1650 | for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld]; |
606 | /* nr_o <- mynorm(A[i]'*m) */ | ||
607 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
150 | PetscCall(SlepcMatDenseMult(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE)); |
608 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
|
150 | nr = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0)); |
609 |
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.
|
150 | PetscCheck(PetscAbs(nr)>PETSC_MACHINE_EPSILON,PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"Linear dependency detected"); |
610 |
6/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 10 times.
|
150 | if (PetscAbs(nr) > 0.7*PetscAbs(nr_o)) break; |
611 | 60 | nr_o = nr; | |
612 | } | ||
613 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
|
100 | ns_[i] = PetscSign(nr); |
614 | /* A[i] <- A[i]/|nr| */ | ||
615 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
100 | nr_abs = PetscAbs(nr); |
616 |
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.
|
100 | PetscCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&nr_abs,&done,&rA_,&one,A+i*ld,&ld_,&info)); |
617 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
100 | SlepcCheckLapackInfo("lascl",info); |
618 | } | ||
619 |
4/6✓ Branch 0 taken 2 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(MatDenseRestoreArray(ds->omat[mat],&A_)); |
620 |
4/6✓ Branch 0 taken 2 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(PetscLogEventEnd(DS_Other,ds,0,0,0)); |
621 |
2/4✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
|
10 | PetscCall(PetscObjectStateIncrease((PetscObject)ds)); |
622 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
10 | if (lindcols) *lindcols = cols; |
623 |
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.
|
2 | PetscFunctionReturn(PETSC_SUCCESS); |
624 | } | ||
625 |