GCC Code Coverage Report


Directory: ./
File: src/sys/classes/bv/tests/test4.c
Date: 2025-10-03 04:28:47
Exec Total Coverage
Lines: 96 109 88.1%
Functions: 1 1 100.0%
Branches: 296 546 54.2%

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 static char help[] = "Test BV operations, changing the number of active columns.\n\n";
12
13 #include <slepcbv.h>
14
15 98 int main(int argc,char **argv)
16 {
17 98 Vec t,v;
18 98 Mat Q,M;
19 98 BV X,Y;
20 98 PetscInt i,j,n=10,kx=6,lx=3,ky=5,ly=2;
21 98 PetscScalar *q,*z;
22 98 PetscReal nrm;
23 98 PetscViewer view;
24 98 PetscBool verbose,trans;
25
26
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
98 PetscFunctionBeginUser;
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.
98 PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
28
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
29
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(PetscOptionsGetInt(NULL,NULL,"-kx",&kx,NULL));
30
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(PetscOptionsGetInt(NULL,NULL,"-lx",&lx,NULL));
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.
98 PetscCall(PetscOptionsGetInt(NULL,NULL,"-ky",&ky,NULL));
32
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(PetscOptionsGetInt(NULL,NULL,"-ly",&ly,NULL));
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.
98 PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
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.
98 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"First BV with %" PetscInt_FMT " active columns (%" PetscInt_FMT " leading columns) of dimension %" PetscInt_FMT ".\n",kx,lx,n));
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.
98 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Second BV with %" PetscInt_FMT " active columns (%" PetscInt_FMT " leading columns) of dimension %" PetscInt_FMT ".\n",ky,ly,n));
36
37 /* Create template vector */
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.
98 PetscCall(VecCreate(PETSC_COMM_WORLD,&t));
39
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(VecSetSizes(t,PETSC_DECIDE,n));
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.
98 PetscCall(VecSetFromOptions(t));
41
42 /* Create BV object X */
43
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(BVCreate(PETSC_COMM_WORLD,&X));
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.
98 PetscCall(PetscObjectSetName((PetscObject)X,"X"));
45
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(BVSetSizesFromVec(X,t,kx+2)); /* two extra columns to test active columns */
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.
98 PetscCall(BVSetFromOptions(X));
47
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(BVSetActiveColumns(X,lx,kx));
48
49 /* Set up viewer */
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.
98 PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view));
51
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.
98 if (verbose) PetscCall(PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB));
52
53 /* Fill X entries */
54
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
1470 for (j=0;j<kx+2;j++) {
55
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1372 PetscCall(BVGetColumn(X,j,&v));
56
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1372 PetscCall(VecSet(v,0.0));
57
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
6860 for (i=0;i<4;i++) {
58
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.
5488 if (i+j<n) PetscCall(VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES));
59 }
60
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1372 PetscCall(VecAssemblyBegin(v));
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.
1372 PetscCall(VecAssemblyEnd(v));
62
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1372 PetscCall(BVRestoreColumn(X,j,&v));
63 }
64
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.
98 if (verbose) PetscCall(BVView(X,view));
65
66 /* Create BV object Y */
67
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(BVCreate(PETSC_COMM_WORLD,&Y));
68
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(PetscObjectSetName((PetscObject)Y,"Y"));
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.
98 PetscCall(BVSetSizesFromVec(Y,t,ky+1));
70
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(BVSetFromOptions(Y));
71
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(BVSetActiveColumns(Y,ly,ky));
72
73 /* Fill Y entries */
74
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
980 for (j=0;j<ky+1;j++) {
75
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
882 PetscCall(BVGetColumn(Y,j,&v));
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.
882 PetscCall(VecSet(v,(PetscScalar)(j+1)/4.0));
77
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
882 PetscCall(BVRestoreColumn(Y,j,&v));
78 }
79
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.
98 if (verbose) PetscCall(BVView(Y,view));
80
81 /* Create Mat */
82
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,kx,ky,NULL,&Q));
83
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(PetscObjectSetName((PetscObject)Q,"Q"));
84
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(MatDenseGetArray(Q,&q));
85
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
1274 for (i=0;i<kx;i++)
86
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
10584 for (j=0;j<ky;j++)
87
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
16072 q[i+j*kx] = (i<j)? 2.0: -0.5;
88
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(MatDenseRestoreArray(Q,&q));
89
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.
98 if (verbose) PetscCall(MatView(Q,NULL));
90
91 /* Test BVResize */
92
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(BVResize(X,kx+4,PETSC_TRUE));
93
94 /* Test BVMult */
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.
98 PetscCall(BVMult(Y,2.0,0.5,X,Q));
96
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
98 if (verbose) {
97 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After BVMult - - - - - - - - -\n"));
98 PetscCall(BVView(Y,view));
99 }
100
101 /* Test BVMultVec */
102
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(BVGetColumn(Y,0,&v));
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.
98 PetscCall(PetscMalloc1(kx-lx,&z));
104 98 z[0] = 2.0;
105
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
882 for (i=1;i<kx-lx;i++) z[i] = -0.5*z[i-1];
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.
98 PetscCall(BVMultVec(X,-1.0,1.0,v,z));
107
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.
98 PetscCall(PetscFree(z));
108
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(BVRestoreColumn(Y,0,&v));
109
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
98 if (verbose) {
110 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After BVMultVec - - - - - - -\n"));
111 PetscCall(BVView(Y,view));
112 }
113
114 /* Test BVDot */
115
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&M));
116
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(PetscObjectSetName((PetscObject)M,"M"));
117
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(BVDot(X,Y,M));
118
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
98 if (verbose) {
119 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After BVDot - - - - - - - - -\n"));
120 PetscCall(MatView(M,NULL));
121 }
122
123 /* Test BVDotVec */
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.
98 PetscCall(BVGetColumn(Y,0,&v));
125
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(PetscMalloc1(kx-lx,&z));
126
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(BVDotVec(X,v,z));
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.
98 PetscCall(BVRestoreColumn(Y,0,&v));
128
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
98 if (verbose) {
129 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After BVDotVec - - - - - - -\n"));
130 PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,kx-lx,z,&v));
131 PetscCall(PetscObjectSetName((PetscObject)v,"z"));
132 PetscCall(VecView(v,view));
133 PetscCall(VecDestroy(&v));
134 }
135
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.
98 PetscCall(PetscFree(z));
136
137 /* Test BVMultInPlace and BVScale */
138
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(PetscOptionsHasName(NULL,NULL,"-trans",&trans));
139
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
98 if (trans) {
140 40 Mat Qt;
141
4/6
✓ Branch 0 taken 2 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(MatTranspose(Q,MAT_INITIAL_MATRIX,&Qt));
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.
40 PetscCall(BVMultInPlaceHermitianTranspose(X,Qt,lx+1,ky));
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.
40 PetscCall(MatDestroy(&Qt));
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.
58 } else PetscCall(BVMultInPlace(X,Q,lx+1,ky));
145
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(BVScale(X,2.0));
146
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
98 if (verbose) {
147 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After BVMultInPlace - - - - -\n"));
148 PetscCall(BVView(X,view));
149 }
150
151 /* Test BVNorm */
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.
98 PetscCall(BVNormColumn(X,lx,NORM_2,&nrm));
153
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"2-Norm of X[%" PetscInt_FMT "] = %g\n",lx,(double)nrm));
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.
98 PetscCall(BVNorm(X,NORM_FROBENIUS,&nrm));
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.
98 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Frobenius Norm of X = %g\n",(double)nrm));
156
157
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(BVDestroy(&X));
158
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(BVDestroy(&Y));
159
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(MatDestroy(&Q));
160
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(MatDestroy(&M));
161
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
98 PetscCall(VecDestroy(&t));
162
3/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
98 PetscCall(SlepcFinalize());
163 return 0;
164 }
165
166 /*TEST
167
168 testset:
169 output_file: output/test4_1.out
170 args: -n 18 -kx 12 -ky 8
171 test:
172 suffix: 1
173 args: -bv_type {{vecs contiguous svec mat}shared output}
174 test:
175 suffix: 1_vecs_vmip
176 args: -bv_type vecs -bv_vecs_vmip 1
177 test:
178 suffix: 1_cuda
179 args: -bv_type {{svec mat}} -vec_type cuda
180 requires: cuda
181 test:
182 suffix: 1_hip
183 args: -bv_type {{svec mat}} -vec_type hip
184 requires: hip
185 test:
186 suffix: 2
187 args: -bv_type {{vecs contiguous svec mat}shared output} -trans
188
189 TEST*/
190