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 | BV implemented as an array of independent Vecs | ||
12 | */ | ||
13 | |||
14 | #include <slepc/private/bvimpl.h> | ||
15 | |||
16 | typedef struct { | ||
17 | Vec *V; | ||
18 | PetscInt vmip; /* Version of BVMultInPlace: | ||
19 | 0: memory-efficient version, uses VecGetArray (default in CPU) | ||
20 | 1: version that allocates (e-s) work vectors in every call (default in GPU) */ | ||
21 | } BV_VECS; | ||
22 | |||
23 | 1008 | static PetscErrorCode BVMult_Vecs(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q) | |
24 | { | ||
25 | 1008 | BV_VECS *y = (BV_VECS*)Y->data,*x = (BV_VECS*)X->data; | |
26 | 1008 | PetscScalar *s=NULL; | |
27 | 1008 | const PetscScalar *q; | |
28 | 1008 | PetscInt i,j,ldq; | |
29 | 1008 | PetscBool trivial=(alpha==1.0)?PETSC_TRUE:PETSC_FALSE; | |
30 | |||
31 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1008 | PetscFunctionBegin; |
32 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1008 | if (Q) { |
33 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
808 | PetscCall(MatDenseGetLDA(Q,&ldq)); |
34 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
808 | if (!trivial) { |
35 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
808 | PetscCall(BVAllocateWork_Private(Y,X->k-X->l)); |
36 | 808 | s = Y->work; | |
37 | } | ||
38 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
808 | PetscCall(MatDenseGetArrayRead(Q,&q)); |
39 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
5372 | for (j=Y->l;j<Y->k;j++) { |
40 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
4564 | PetscCall(VecScale(y->V[Y->nc+j],beta)); |
41 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
4564 | if (!trivial) { |
42 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
27344 | for (i=X->l;i<X->k;i++) s[i-X->l] = alpha*q[i+j*ldq]; |
43 | ✗ | } else s = (PetscScalar*)(q+j*ldq+X->l); | |
44 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
4564 | PetscCall(VecMAXPY(y->V[Y->nc+j],X->k-X->l,s,x->V+X->nc+X->l)); |
45 | } | ||
46 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
808 | PetscCall(MatDenseRestoreArrayRead(Q,&q)); |
47 | } else { | ||
48 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1000 | for (j=0;j<Y->k-Y->l;j++) { |
49 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
800 | PetscCall(VecScale(y->V[Y->nc+Y->l+j],beta)); |
50 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
800 | PetscCall(VecAXPY(y->V[Y->nc+Y->l+j],alpha,x->V[X->nc+X->l+j])); |
51 | } | ||
52 | } | ||
53 |
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.
|
198 | PetscFunctionReturn(PETSC_SUCCESS); |
54 | } | ||
55 | |||
56 | 9674 | static PetscErrorCode BVMultVec_Vecs(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q) | |
57 | { | ||
58 | 9674 | BV_VECS *x = (BV_VECS*)X->data; | |
59 | 9674 | PetscScalar *s=NULL,*qq=q; | |
60 | 9674 | PetscInt i; | |
61 | 9674 | PetscBool trivial=(alpha==1.0)?PETSC_TRUE:PETSC_FALSE; | |
62 | |||
63 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
9674 | PetscFunctionBegin; |
64 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
9674 | if (!trivial) { |
65 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
9674 | PetscCall(BVAllocateWork_Private(X,X->k-X->l)); |
66 | 9674 | s = X->work; | |
67 | } | ||
68 |
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.
|
9674 | if (!q) PetscCall(VecGetArray(X->buffer,&qq)); |
69 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
9674 | PetscCall(VecScale(y,beta)); |
70 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
9674 | if (!trivial) { |
71 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
72978 | for (i=0;i<X->k-X->l;i++) s[i] = alpha*qq[i]; |
72 | ✗ | } else s = qq; | |
73 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
9674 | PetscCall(VecMAXPY(y,X->k-X->l,s,x->V+X->nc+X->l)); |
74 |
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.
|
9674 | if (!q) PetscCall(VecRestoreArray(X->buffer,&qq)); |
75 |
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.
|
1883 | PetscFunctionReturn(PETSC_SUCCESS); |
76 | } | ||
77 | |||
78 | /* | ||
79 | BVMultInPlace_Vecs_ME - V(:,s:e-1) = V*Q(:,s:e-1) for regular vectors. | ||
80 | |||
81 | Memory-efficient version, uses VecGetArray (default in CPU) | ||
82 | |||
83 | Writing V = [ V1 V2 V3 ] and Q(:,s:e-1) = [ Q1 Q2 Q3 ]', where V2 | ||
84 | corresponds to the columns s:e-1, the computation is done as | ||
85 | V2 := V2*Q2 + V1*Q1 + V3*Q3 | ||
86 | */ | ||
87 | 1160 | static PetscErrorCode BVMultInPlace_Vecs_ME(BV V,Mat Q,PetscInt s,PetscInt e) | |
88 | { | ||
89 | 1160 | BV_VECS *ctx = (BV_VECS*)V->data; | |
90 | 1160 | const PetscScalar *q; | |
91 | 1160 | PetscInt i,ldq; | |
92 | |||
93 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1160 | PetscFunctionBegin; |
94 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1160 | PetscCall(MatDenseGetLDA(Q,&ldq)); |
95 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1160 | PetscCall(MatDenseGetArrayRead(Q,&q)); |
96 | /* V2 := V2*Q2 */ | ||
97 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1160 | PetscCall(BVMultInPlace_Vecs_Private(V,V->n,e-s,ldq,ctx->V+V->nc+s,q+s*ldq+s,PETSC_FALSE)); |
98 | /* V2 += V1*Q1 + V3*Q3 */ | ||
99 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
5642 | for (i=s;i<e;i++) { |
100 |
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.
|
4482 | if (PetscUnlikely(s>V->l)) PetscCall(VecMAXPY(ctx->V[V->nc+i],s-V->l,q+i*ldq+V->l,ctx->V+V->nc+V->l)); |
101 |
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.
|
4482 | if (V->k>e) PetscCall(VecMAXPY(ctx->V[V->nc+i],V->k-e,q+i*ldq+e,ctx->V+V->nc+e)); |
102 | } | ||
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.
|
1160 | PetscCall(MatDenseRestoreArrayRead(Q,&q)); |
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.
|
224 | PetscFunctionReturn(PETSC_SUCCESS); |
105 | } | ||
106 | |||
107 | /* | ||
108 | BVMultInPlace_Vecs_Alloc - V(:,s:e-1) = V*Q(:,s:e-1) for regular vectors. | ||
109 | |||
110 | Version that allocates (e-s) work vectors in every call (default in GPU) | ||
111 | */ | ||
112 | 18 | static PetscErrorCode BVMultInPlace_Vecs_Alloc(BV V,Mat Q,PetscInt s,PetscInt e) | |
113 | { | ||
114 | 18 | BV_VECS *ctx = (BV_VECS*)V->data; | |
115 | 18 | const PetscScalar *q; | |
116 | 18 | PetscInt i,ldq; | |
117 | 18 | Vec *W,z; | |
118 | |||
119 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
18 | PetscFunctionBegin; |
120 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
18 | PetscCall(MatDenseGetLDA(Q,&ldq)); |
121 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
18 | PetscCall(MatDenseGetArrayRead(Q,&q)); |
122 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
18 | PetscCall(BVCreateVec(V,&z)); |
123 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
18 | PetscCall(VecDuplicateVecs(z,e-s,&W)); |
124 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
18 | PetscCall(VecDestroy(&z)); |
125 |
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.
|
98 | for (i=s;i<e;i++) PetscCall(VecMAXPY(W[i-s],V->k-V->l,q+i*ldq+V->l,ctx->V+V->nc+V->l)); |
126 |
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.
|
98 | for (i=s;i<e;i++) PetscCall(VecCopy(W[i-s],ctx->V[V->nc+i])); |
127 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
18 | PetscCall(VecDestroyVecs(e-s,&W)); |
128 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
18 | PetscCall(MatDenseRestoreArrayRead(Q,&q)); |
129 |
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); |
130 | } | ||
131 | |||
132 | /* | ||
133 | BVMultInPlaceHermitianTranspose_Vecs - V(:,s:e-1) = V*Q'(:,s:e-1) for regular vectors. | ||
134 | */ | ||
135 | 10 | static PetscErrorCode BVMultInPlaceHermitianTranspose_Vecs(BV V,Mat Q,PetscInt s,PetscInt e) | |
136 | { | ||
137 | 10 | BV_VECS *ctx = (BV_VECS*)V->data; | |
138 | 10 | const PetscScalar *q; | |
139 | 10 | PetscInt i,j,ldq,n; | |
140 | |||
141 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
142 |
4/6✓ Branch 0 taken 2 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(MatGetSize(Q,NULL,&n)); |
143 |
4/6✓ Branch 0 taken 2 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(MatDenseGetLDA(Q,&ldq)); |
144 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatDenseGetArrayRead(Q,&q)); |
145 | /* V2 := V2*Q2' */ | ||
146 |
4/6✓ Branch 0 taken 2 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(BVMultInPlace_Vecs_Private(V,V->n,e-s,ldq,ctx->V+V->nc+s,q+s*ldq+s,PETSC_TRUE)); |
147 | /* V2 += V1*Q1' + V3*Q3' */ | ||
148 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
50 | for (i=s;i<e;i++) { |
149 |
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.
|
80 | for (j=V->l;j<s;j++) PetscCall(VecAXPY(ctx->V[V->nc+i],q[i+j*ldq],ctx->V[V->nc+j])); |
150 |
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.
|
200 | for (j=e;j<n;j++) PetscCall(VecAXPY(ctx->V[V->nc+i],q[i+j*ldq],ctx->V[V->nc+j])); |
151 | } | ||
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.
|
10 | PetscCall(MatDenseRestoreArrayRead(Q,&q)); |
153 |
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); |
154 | } | ||
155 | |||
156 | 1558 | static PetscErrorCode BVDot_Vecs(BV X,BV Y,Mat M) | |
157 | { | ||
158 | 1558 | BV_VECS *x = (BV_VECS*)X->data,*y = (BV_VECS*)Y->data; | |
159 | 1558 | PetscScalar *m; | |
160 | 1558 | PetscInt j,ldm; | |
161 | |||
162 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1558 | PetscFunctionBegin; |
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.
|
1558 | PetscCall(MatDenseGetLDA(M,&ldm)); |
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.
|
1558 | PetscCall(MatDenseGetArray(M,&m)); |
165 |
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.
|
10378 | for (j=X->l;j<X->k;j++) PetscCall(VecMDot(x->V[X->nc+j],Y->k-Y->l,y->V+Y->nc+Y->l,m+j*ldm+Y->l)); |
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.
|
1558 | PetscCall(MatDenseRestoreArray(M,&m)); |
167 |
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.
|
308 | PetscFunctionReturn(PETSC_SUCCESS); |
168 | } | ||
169 | |||
170 | 13584 | static PetscErrorCode BVDotVec_Vecs(BV X,Vec y,PetscScalar *q) | |
171 | { | ||
172 | 13584 | BV_VECS *x = (BV_VECS*)X->data; | |
173 | 13584 | Vec z = y; | |
174 | 13584 | PetscScalar *qq=q; | |
175 | |||
176 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
13584 | PetscFunctionBegin; |
177 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
13584 | if (PetscUnlikely(X->matrix)) { |
178 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
570 | PetscCall(BV_IPMatMult(X,y)); |
179 | 570 | z = X->Bx; | |
180 | } | ||
181 |
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.
|
13584 | if (!q) PetscCall(VecGetArray(X->buffer,&qq)); |
182 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
13584 | PetscCall(VecMDot(z,X->k-X->l,x->V+X->nc+X->l,qq)); |
183 |
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.
|
13584 | if (!q) PetscCall(VecRestoreArray(X->buffer,&qq)); |
184 |
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.
|
2665 | PetscFunctionReturn(PETSC_SUCCESS); |
185 | } | ||
186 | |||
187 | 40 | static PetscErrorCode BVDotVec_Begin_Vecs(BV X,Vec y,PetscScalar *m) | |
188 | { | ||
189 | 40 | BV_VECS *x = (BV_VECS*)X->data; | |
190 | 40 | Vec z = y; | |
191 | |||
192 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
40 | PetscFunctionBegin; |
193 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
40 | if (PetscUnlikely(X->matrix)) { |
194 | ✗ | PetscCall(BV_IPMatMult(X,y)); | |
195 | ✗ | z = X->Bx; | |
196 | } | ||
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.
|
40 | PetscCall(VecMDotBegin(z,X->k-X->l,x->V+X->nc+X->l,m)); |
198 |
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.
|
8 | PetscFunctionReturn(PETSC_SUCCESS); |
199 | } | ||
200 | |||
201 | 40 | static PetscErrorCode BVDotVec_End_Vecs(BV X,Vec y,PetscScalar *m) | |
202 | { | ||
203 | 40 | BV_VECS *x = (BV_VECS*)X->data; | |
204 | |||
205 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
40 | PetscFunctionBegin; |
206 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(VecMDotEnd(y,X->k-X->l,x->V+X->nc+X->l,m)); |
207 |
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.
|
8 | PetscFunctionReturn(PETSC_SUCCESS); |
208 | } | ||
209 | |||
210 | 6315 | static PetscErrorCode BVScale_Vecs(BV bv,PetscInt j,PetscScalar alpha) | |
211 | { | ||
212 | 6315 | PetscInt i; | |
213 | 6315 | BV_VECS *ctx = (BV_VECS*)bv->data; | |
214 | |||
215 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
6315 | PetscFunctionBegin; |
216 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
6315 | if (PetscUnlikely(j<0)) { |
217 |
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.
|
528 | for (i=bv->l;i<bv->k;i++) PetscCall(VecScale(ctx->V[bv->nc+i],alpha)); |
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.
|
6247 | } else PetscCall(VecScale(ctx->V[bv->nc+j],alpha)); |
219 |
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.
|
1229 | PetscFunctionReturn(PETSC_SUCCESS); |
220 | } | ||
221 | |||
222 | 1156 | static PetscErrorCode BVNorm_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val) | |
223 | { | ||
224 | 1156 | PetscInt i; | |
225 | 1156 | PetscReal nrm; | |
226 | 1156 | BV_VECS *ctx = (BV_VECS*)bv->data; | |
227 | |||
228 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1156 | PetscFunctionBegin; |
229 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1156 | if (PetscUnlikely(j<0)) { |
230 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
708 | PetscCheck(type==NORM_FROBENIUS,PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not implemented in BVVECS"); |
231 | 708 | *val = 0.0; | |
232 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
4428 | for (i=bv->l;i<bv->k;i++) { |
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.
|
3720 | PetscCall(VecNorm(ctx->V[bv->nc+i],NORM_2,&nrm)); |
234 | 3720 | *val += nrm*nrm; | |
235 | } | ||
236 | 708 | *val = PetscSqrtReal(*val); | |
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.
|
448 | } else PetscCall(VecNorm(ctx->V[bv->nc+j],type,val)); |
238 |
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.
|
222 | PetscFunctionReturn(PETSC_SUCCESS); |
239 | } | ||
240 | |||
241 | 20 | static PetscErrorCode BVNorm_Begin_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val) | |
242 | { | ||
243 | 20 | BV_VECS *ctx = (BV_VECS*)bv->data; | |
244 | |||
245 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
20 | PetscFunctionBegin; |
246 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
20 | PetscCheck(j>=0,PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not implemented in BVVECS"); |
247 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(VecNormBegin(ctx->V[bv->nc+j],type,val)); |
248 |
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.
|
4 | PetscFunctionReturn(PETSC_SUCCESS); |
249 | } | ||
250 | |||
251 | 20 | static PetscErrorCode BVNorm_End_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val) | |
252 | { | ||
253 | 20 | BV_VECS *ctx = (BV_VECS*)bv->data; | |
254 | |||
255 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
20 | PetscFunctionBegin; |
256 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
20 | PetscCheck(j>=0,PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not implemented in BVVECS"); |
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.
|
20 | PetscCall(VecNormEnd(ctx->V[bv->nc+j],type,val)); |
258 |
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.
|
4 | PetscFunctionReturn(PETSC_SUCCESS); |
259 | } | ||
260 | |||
261 | 45 | static PetscErrorCode BVNormalize_Vecs(BV bv,PetscScalar *eigi) | |
262 | { | ||
263 | 45 | BV_VECS *ctx = (BV_VECS*)bv->data; | |
264 | 45 | PetscInt i; | |
265 | |||
266 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
45 | PetscFunctionBegin; |
267 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
420 | for (i=bv->l;i<bv->k;i++) { |
268 | #if !defined(PETSC_USE_COMPLEX) | ||
269 |
4/4✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
|
240 | if (eigi && eigi[i] != 0.0) { |
270 |
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.
|
30 | PetscCall(VecNormalizeComplex(ctx->V[bv->nc+i],ctx->V[bv->nc+i+1],PETSC_TRUE,NULL)); |
271 | 30 | i++; | |
272 | } else | ||
273 | #endif | ||
274 | { | ||
275 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
375 | PetscCall(VecNormalize(ctx->V[bv->nc+i],NULL)); |
276 | } | ||
277 | } | ||
278 |
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.
|
9 | PetscFunctionReturn(PETSC_SUCCESS); |
279 | } | ||
280 | |||
281 | 940 | static PetscErrorCode BVMatMult_Vecs(BV V,Mat A,BV W) | |
282 | { | ||
283 | 940 | BV_VECS *v = (BV_VECS*)V->data,*w = (BV_VECS*)W->data; | |
284 | 940 | PetscInt j; | |
285 | 940 | Mat Vmat,Wmat; | |
286 | |||
287 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
940 | PetscFunctionBegin; |
288 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
940 | if (V->vmm) { |
289 |
4/6✓ Branch 0 taken 2 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 | PetscCall(BVGetMat(V,&Vmat)); |
290 |
4/6✓ Branch 0 taken 2 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 | PetscCall(BVGetMat(W,&Wmat)); |
291 |
4/6✓ Branch 0 taken 2 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 | PetscCall(MatProductCreateWithMat(A,Vmat,NULL,Wmat)); |
292 |
4/6✓ Branch 0 taken 2 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 | PetscCall(MatProductSetType(Wmat,MATPRODUCT_AB)); |
293 |
4/6✓ Branch 0 taken 2 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 | PetscCall(MatProductSetFromOptions(Wmat)); |
294 |
4/6✓ Branch 0 taken 2 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 | PetscCall(MatProductSymbolic(Wmat)); |
295 |
4/6✓ Branch 0 taken 2 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 | PetscCall(MatProductNumeric(Wmat)); |
296 |
4/6✓ Branch 0 taken 2 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 | PetscCall(MatProductClear(Wmat)); |
297 |
4/6✓ Branch 0 taken 2 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 | PetscCall(BVRestoreMat(V,&Vmat)); |
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.
|
50 | PetscCall(BVRestoreMat(W,&Wmat)); |
299 | } else { | ||
300 |
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.
|
5130 | for (j=0;j<V->k-V->l;j++) PetscCall(MatMult(A,v->V[V->nc+V->l+j],w->V[W->nc+W->l+j])); |
301 | } | ||
302 |
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.
|
188 | PetscFunctionReturn(PETSC_SUCCESS); |
303 | } | ||
304 | |||
305 | 970 | static PetscErrorCode BVCopy_Vecs(BV V,BV W) | |
306 | { | ||
307 | 970 | BV_VECS *v = (BV_VECS*)V->data,*w = (BV_VECS*)W->data; | |
308 | 970 | PetscInt j; | |
309 | |||
310 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
970 | PetscFunctionBegin; |
311 |
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.
|
7260 | for (j=0;j<V->k-V->l;j++) PetscCall(VecCopy(v->V[V->nc+V->l+j],w->V[W->nc+W->l+j])); |
312 |
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.
|
194 | PetscFunctionReturn(PETSC_SUCCESS); |
313 | } | ||
314 | |||
315 | 692 | static PetscErrorCode BVCopyColumn_Vecs(BV V,PetscInt j,PetscInt i) | |
316 | { | ||
317 | 692 | BV_VECS *v = (BV_VECS*)V->data; | |
318 | |||
319 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
692 | PetscFunctionBegin; |
320 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
692 | PetscCall(VecCopy(v->V[V->nc+j],v->V[V->nc+i])); |
321 |
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.
|
134 | PetscFunctionReturn(PETSC_SUCCESS); |
322 | } | ||
323 | |||
324 | 90 | static PetscErrorCode BVResize_Vecs(BV bv,PetscInt m,PetscBool copy) | |
325 | { | ||
326 | 90 | BV_VECS *ctx = (BV_VECS*)bv->data; | |
327 | 90 | Vec *newV,z; | |
328 | 90 | PetscInt j; | |
329 | 90 | char str[50]; | |
330 | |||
331 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
90 | PetscFunctionBegin; |
332 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
90 | PetscCall(BVCreateVec(bv,&z)); |
333 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
90 | PetscCall(VecDuplicateVecs(z,m,&newV)); |
334 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
90 | PetscCall(VecDestroy(&z)); |
335 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
90 | if (((PetscObject)bv)->name) { |
336 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1110 | for (j=0;j<m;j++) { |
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.
|
1020 | PetscCall(PetscSNPrintf(str,sizeof(str),"%s_%" PetscInt_FMT,((PetscObject)bv)->name,j)); |
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.
|
1020 | PetscCall(PetscObjectSetName((PetscObject)newV[j],str)); |
339 | } | ||
340 | } | ||
341 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
90 | if (copy) { |
342 |
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.
|
450 | for (j=0;j<PetscMin(m,bv->m);j++) PetscCall(VecCopy(ctx->V[j],newV[j])); |
343 | } | ||
344 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
90 | PetscCall(VecDestroyVecs(bv->m,&ctx->V)); |
345 | 90 | ctx->V = newV; | |
346 |
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.
|
90 | PetscFunctionReturn(PETSC_SUCCESS); |
347 | } | ||
348 | |||
349 | 46207 | static PetscErrorCode BVGetColumn_Vecs(BV bv,PetscInt j,Vec *v) | |
350 | { | ||
351 | 46207 | BV_VECS *ctx = (BV_VECS*)bv->data; | |
352 | 46207 | PetscInt l; | |
353 | |||
354 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
46207 | PetscFunctionBegin; |
355 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
|
46207 | l = BVAvailableVec; |
356 | 46207 | bv->cv[l] = ctx->V[bv->nc+j]; | |
357 |
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.
|
46207 | PetscFunctionReturn(PETSC_SUCCESS); |
358 | } | ||
359 | |||
360 | 46207 | static PetscErrorCode BVRestoreColumn_Vecs(BV bv,PetscInt j,Vec *v) | |
361 | { | ||
362 | 46207 | PetscInt l; | |
363 | |||
364 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
46207 | PetscFunctionBegin; |
365 | 46207 | l = (j==bv->ci[0])? 0: 1; | |
366 | 46207 | bv->cv[l] = NULL; | |
367 |
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.
|
46207 | PetscFunctionReturn(PETSC_SUCCESS); |
368 | } | ||
369 | |||
370 | 590 | static PetscErrorCode BVGetArray_Vecs(BV bv,PetscScalar **a) | |
371 | { | ||
372 | 590 | BV_VECS *ctx = (BV_VECS*)bv->data; | |
373 | 590 | PetscInt j; | |
374 | 590 | const PetscScalar *p; | |
375 | |||
376 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
590 | PetscFunctionBegin; |
377 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
590 | PetscCall(PetscMalloc1((bv->nc+bv->m)*bv->n,a)); |
378 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
5480 | for (j=0;j<bv->nc+bv->m;j++) { |
379 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
4890 | PetscCall(VecGetArrayRead(ctx->V[j],&p)); |
380 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
4890 | PetscCall(PetscArraycpy(*a+j*bv->n,p,bv->n)); |
381 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
4890 | PetscCall(VecRestoreArrayRead(ctx->V[j],&p)); |
382 | } | ||
383 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
118 | PetscFunctionReturn(PETSC_SUCCESS); |
384 | } | ||
385 | |||
386 | 590 | static PetscErrorCode BVRestoreArray_Vecs(BV bv,PetscScalar **a) | |
387 | { | ||
388 | 590 | BV_VECS *ctx = (BV_VECS*)bv->data; | |
389 | 590 | PetscInt j; | |
390 | 590 | PetscScalar *p; | |
391 | |||
392 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
590 | PetscFunctionBegin; |
393 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
5480 | for (j=0;j<bv->nc+bv->m;j++) { |
394 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
4890 | PetscCall(VecGetArray(ctx->V[j],&p)); |
395 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
4890 | PetscCall(PetscArraycpy(p,*a+j*bv->n,bv->n)); |
396 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
4890 | PetscCall(VecRestoreArray(ctx->V[j],&p)); |
397 | } | ||
398 |
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.
|
590 | PetscCall(PetscFree(*a)); |
399 |
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.
|
118 | PetscFunctionReturn(PETSC_SUCCESS); |
400 | } | ||
401 | |||
402 | 80 | static PetscErrorCode BVGetArrayRead_Vecs(BV bv,const PetscScalar **a) | |
403 | { | ||
404 | 80 | BV_VECS *ctx = (BV_VECS*)bv->data; | |
405 | 80 | PetscInt j; | |
406 | 80 | const PetscScalar *p; | |
407 | |||
408 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
80 | PetscFunctionBegin; |
409 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
80 | PetscCall(PetscMalloc1((bv->nc+bv->m)*bv->n,(PetscScalar**)a)); |
410 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
540 | for (j=0;j<bv->nc+bv->m;j++) { |
411 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
460 | PetscCall(VecGetArrayRead(ctx->V[j],&p)); |
412 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
460 | PetscCall(PetscArraycpy((PetscScalar*)*a+j*bv->n,p,bv->n)); |
413 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
460 | PetscCall(VecRestoreArrayRead(ctx->V[j],&p)); |
414 | } | ||
415 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
16 | PetscFunctionReturn(PETSC_SUCCESS); |
416 | } | ||
417 | |||
418 | 80 | static PetscErrorCode BVRestoreArrayRead_Vecs(BV bv,const PetscScalar **a) | |
419 | { | ||
420 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
80 | PetscFunctionBegin; |
421 |
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.
|
80 | PetscCall(PetscFree(*a)); |
422 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
16 | PetscFunctionReturn(PETSC_SUCCESS); |
423 | } | ||
424 | |||
425 | /* | ||
426 | Sets the value of vmip flag and resets ops->multinplace accordingly | ||
427 | */ | ||
428 | 6488 | static inline PetscErrorCode BVVecsSetVmip(BV bv,PetscInt vmip) | |
429 | { | ||
430 | 6488 | typedef PetscErrorCode (*fmultinplace)(BV,Mat,PetscInt,PetscInt); | |
431 | 6488 | fmultinplace multinplace[2] = {BVMultInPlace_Vecs_ME, BVMultInPlace_Vecs_Alloc}; | |
432 | 6488 | BV_VECS *ctx = (BV_VECS*)bv->data; | |
433 | |||
434 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
6488 | PetscFunctionBegin; |
435 | 6488 | ctx->vmip = vmip; | |
436 | 6488 | bv->ops->multinplace = multinplace[vmip]; | |
437 |
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.
|
6488 | PetscFunctionReturn(PETSC_SUCCESS); |
438 | } | ||
439 | |||
440 | 1174 | static PetscErrorCode BVSetFromOptions_Vecs(BV bv,PetscOptionItems PetscOptionsObject) | |
441 | { | ||
442 | 1174 | BV_VECS *ctx = (BV_VECS*)bv->data; | |
443 | |||
444 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1174 | PetscFunctionBegin; |
445 |
1/12✗ 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.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
1174 | PetscOptionsHeadBegin(PetscOptionsObject,"BV Vecs Options"); |
446 | |||
447 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1174 | PetscCall(PetscOptionsRangeInt("-bv_vecs_vmip","Version of BVMultInPlace operation","",ctx->vmip,&ctx->vmip,NULL,0,1)); |
448 |
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.
|
1174 | PetscCall(BVVecsSetVmip(bv,ctx->vmip)); |
449 | |||
450 |
1/14✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
1174 | PetscOptionsHeadEnd(); |
451 |
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.
|
226 | PetscFunctionReturn(PETSC_SUCCESS); |
452 | } | ||
453 | |||
454 | 120 | PetscErrorCode BVView_Vecs(BV bv,PetscViewer viewer) | |
455 | { | ||
456 | 120 | BV_VECS *ctx = (BV_VECS*)bv->data; | |
457 | 120 | PetscInt j; | |
458 | 120 | PetscViewerFormat format; | |
459 | 120 | PetscBool isascii,ismatlab=PETSC_FALSE; | |
460 | 120 | const char *bvname,*name; | |
461 | |||
462 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
120 | PetscFunctionBegin; |
463 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); |
464 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
120 | if (isascii) { |
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.
|
120 | PetscCall(PetscViewerGetFormat(viewer,&format)); |
466 |
8/14✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 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.
|
120 | if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(PETSC_SUCCESS); |
467 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
100 | if (format == PETSC_VIEWER_ASCII_MATLAB) ismatlab = PETSC_TRUE; |
468 | } | ||
469 | ✗ | if (ismatlab) { | |
470 | ✗ | PetscCall(PetscObjectGetName((PetscObject)bv,&bvname)); | |
471 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer,"%s=[];\n",bvname)); | |
472 | } | ||
473 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
480 | for (j=bv->nc;j<bv->nc+bv->m;j++) { |
474 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
380 | PetscCall(VecView(ctx->V[j],viewer)); |
475 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
380 | if (ismatlab) { |
476 | ✗ | PetscCall(PetscObjectGetName((PetscObject)ctx->V[j],&name)); | |
477 |
0/6✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
380 | PetscCall(PetscViewerASCIIPrintf(viewer,"%s=[%s,%s];clear %s\n",bvname,bvname,name,name)); |
478 | } | ||
479 | } | ||
480 |
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.
|
20 | PetscFunctionReturn(PETSC_SUCCESS); |
481 | } | ||
482 | |||
483 | 3244 | static PetscErrorCode BVDestroy_Vecs(BV bv) | |
484 | { | ||
485 | 3244 | BV_VECS *ctx = (BV_VECS*)bv->data; | |
486 | |||
487 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
3244 | PetscFunctionBegin; |
488 |
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.
|
3244 | if (!bv->issplit) PetscCall(VecDestroyVecs(bv->nc+bv->m,&ctx->V)); |
489 |
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.
|
3244 | PetscCall(PetscFree(bv->data)); |
490 |
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.
|
640 | PetscFunctionReturn(PETSC_SUCCESS); |
491 | } | ||
492 | |||
493 | 2070 | static PetscErrorCode BVDuplicate_Vecs(BV V,BV W) | |
494 | { | ||
495 | 2070 | BV_VECS *ctx = (BV_VECS*)V->data; | |
496 | |||
497 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
2070 | PetscFunctionBegin; |
498 |
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.
|
2070 | PetscCall(BVVecsSetVmip(W,ctx->vmip)); |
499 |
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.
|
2070 | PetscFunctionReturn(PETSC_SUCCESS); |
500 | } | ||
501 | |||
502 | 3244 | SLEPC_EXTERN PetscErrorCode BVCreate_Vecs(BV bv) | |
503 | { | ||
504 | 3244 | BV_VECS *ctx; | |
505 | 3244 | PetscInt j,lsplit; | |
506 | 3244 | PetscBool isgpu; | |
507 | 3244 | char str[50]; | |
508 | 3244 | BV parent; | |
509 | 3244 | Vec *Vpar,z; | |
510 | |||
511 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
3244 | PetscFunctionBegin; |
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.
|
3244 | PetscCall(PetscNew(&ctx)); |
513 | 3244 | bv->data = (void*)ctx; | |
514 | |||
515 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
3244 | if (PetscUnlikely(bv->issplit)) { |
516 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
800 | PetscCheck(bv->issplit>0,PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"BVVECS does not support BVGetSplitRows()"); |
517 | /* split BV: share the Vecs of the parent BV */ | ||
518 | 800 | parent = bv->splitparent; | |
519 | 800 | lsplit = parent->lsplit; | |
520 | 800 | Vpar = ((BV_VECS*)parent->data)->V; | |
521 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
800 | ctx->V = (bv->issplit==1)? Vpar: Vpar+lsplit; |
522 | } else { | ||
523 | /* regular BV: create array of Vecs to store the BV columns */ | ||
524 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2444 | PetscCall(BVCreateVec(bv,&z)); |
525 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2444 | PetscCall(VecDuplicateVecs(z,bv->m,&ctx->V)); |
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.
|
2444 | PetscCall(VecDestroy(&z)); |
527 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2444 | if (((PetscObject)bv)->name) { |
528 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
9822 | for (j=0;j<bv->m;j++) { |
529 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8714 | PetscCall(PetscSNPrintf(str,sizeof(str),"%s_%" PetscInt_FMT,((PetscObject)bv)->name,j)); |
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.
|
8714 | PetscCall(PetscObjectSetName((PetscObject)ctx->V[j],str)); |
531 | } | ||
532 | } | ||
533 | } | ||
534 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
3244 | if (!bv->ld) bv->ld = bv->n; |
535 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
3244 | PetscCheck(bv->ld==bv->n,PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"BVVECS does not support a user-defined leading dimension"); |
536 | |||
537 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
3244 | if (PetscUnlikely(bv->Acreate)) { |
538 |
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.
|
180 | for (j=0;j<bv->m;j++) PetscCall(MatGetColumnVector(bv->Acreate,ctx->V[j],j)); |
539 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(MatDestroy(&bv->Acreate)); |
540 | } | ||
541 | |||
542 | /* Default version of BVMultInPlace */ | ||
543 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3244 | PetscCall(PetscStrcmpAny(bv->vtype,&isgpu,VECSEQCUDA,VECMPICUDA,VECSEQHIP,VECMPIHIP,"")); |
544 | 3244 | ctx->vmip = isgpu? 1: 0; | |
545 | |||
546 | /* Default BVMatMult method */ | ||
547 | 3244 | bv->vmm = BV_MATMULT_VECS; | |
548 | |||
549 | /* Deferred call to setfromoptions */ | ||
550 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
3244 | if (bv->defersfo) { |
551 |
8/10✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ 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.
|
30 | PetscObjectOptionsBegin((PetscObject)bv); |
552 |
4/6✓ Branch 0 taken 2 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(BVSetFromOptions_Vecs(bv,PetscOptionsObject)); |
553 |
4/6✓ Branch 0 taken 2 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 | PetscOptionsEnd(); |
554 | } | ||
555 |
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.
|
3244 | PetscCall(BVVecsSetVmip(bv,ctx->vmip)); |
556 | |||
557 | 3244 | bv->ops->mult = BVMult_Vecs; | |
558 | 3244 | bv->ops->multvec = BVMultVec_Vecs; | |
559 | 3244 | bv->ops->multinplacetrans = BVMultInPlaceHermitianTranspose_Vecs; | |
560 | 3244 | bv->ops->dot = BVDot_Vecs; | |
561 | 3244 | bv->ops->dotvec = BVDotVec_Vecs; | |
562 | 3244 | bv->ops->dotvec_begin = BVDotVec_Begin_Vecs; | |
563 | 3244 | bv->ops->dotvec_end = BVDotVec_End_Vecs; | |
564 | 3244 | bv->ops->scale = BVScale_Vecs; | |
565 | 3244 | bv->ops->norm = BVNorm_Vecs; | |
566 | 3244 | bv->ops->norm_begin = BVNorm_Begin_Vecs; | |
567 | 3244 | bv->ops->norm_end = BVNorm_End_Vecs; | |
568 | 3244 | bv->ops->normalize = BVNormalize_Vecs; | |
569 | 3244 | bv->ops->matmult = BVMatMult_Vecs; | |
570 | 3244 | bv->ops->copy = BVCopy_Vecs; | |
571 | 3244 | bv->ops->copycolumn = BVCopyColumn_Vecs; | |
572 | 3244 | bv->ops->resize = BVResize_Vecs; | |
573 | 3244 | bv->ops->getcolumn = BVGetColumn_Vecs; | |
574 | 3244 | bv->ops->restorecolumn = BVRestoreColumn_Vecs; | |
575 | 3244 | bv->ops->getarray = BVGetArray_Vecs; | |
576 | 3244 | bv->ops->restorearray = BVRestoreArray_Vecs; | |
577 | 3244 | bv->ops->getarrayread = BVGetArrayRead_Vecs; | |
578 | 3244 | bv->ops->restorearrayread = BVRestoreArrayRead_Vecs; | |
579 | 3244 | bv->ops->getmat = BVGetMat_Default; | |
580 | 3244 | bv->ops->restoremat = BVRestoreMat_Default; | |
581 | 3244 | bv->ops->destroy = BVDestroy_Vecs; | |
582 | 3244 | bv->ops->duplicate = BVDuplicate_Vecs; | |
583 | 3244 | bv->ops->setfromoptions = BVSetFromOptions_Vecs; | |
584 | 3244 | bv->ops->view = BVView_Vecs; | |
585 |
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.
|
3244 | PetscFunctionReturn(PETSC_SUCCESS); |
586 | } | ||
587 |