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 BVGetSplit().\n\n"; | ||
12 | |||
13 | #include <slepcbv.h> | ||
14 | |||
15 | /* | ||
16 | Print the first row of a BV | ||
17 | */ | ||
18 | 576 | PetscErrorCode PrintFirstRow(BV X) | |
19 | { | ||
20 | 576 | PetscMPIInt rank; | |
21 | 576 | PetscInt i,ld,k,nc; | |
22 | 576 | const PetscScalar *pX; | |
23 | 576 | const char *name; | |
24 | |||
25 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
576 | PetscFunctionBeginUser; |
26 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
576 | PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)X),&rank)); |
27 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
576 | if (!rank) { |
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.
|
288 | PetscCall(BVGetActiveColumns(X,NULL,&k)); |
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.
|
288 | PetscCall(BVGetLeadingDimension(X,&ld)); |
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.
|
288 | PetscCall(BVGetNumConstraints(X,&nc)); |
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.
|
288 | PetscCall(PetscObjectGetName((PetscObject)X,&name)); |
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.
|
288 | PetscCall(PetscPrintf(PetscObjectComm((PetscObject)X),"First row of %s =\n",name)); |
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.
|
288 | PetscCall(BVGetArrayRead(X,&pX)); |
34 |
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.
|
2016 | for (i=0;i<nc+k;i++) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)X),"%g ",(double)PetscRealPart(pX[i*ld]))); |
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.
|
288 | PetscCall(PetscPrintf(PetscObjectComm((PetscObject)X),"\n")); |
36 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
288 | PetscCall(BVRestoreArrayRead(X,&pX)); |
37 | } | ||
38 |
5/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 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
96 | PetscFunctionReturn(PETSC_SUCCESS); |
39 | } | ||
40 | |||
41 | 192 | int main(int argc,char **argv) | |
42 | { | ||
43 | 192 | Vec t,v,*C; | |
44 | 192 | BV X,L,R; | |
45 | 192 | PetscInt i,j,n=10,k=5,l=3,nc=0; | |
46 | 192 | PetscReal norm; | |
47 | 192 | PetscScalar alpha; | |
48 | 192 | PetscViewer view; | |
49 | 192 | PetscBool verbose; | |
50 | |||
51 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
192 | PetscFunctionBeginUser; |
52 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
53 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
54 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL)); |
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.
|
192 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL)); |
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.
|
192 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-nc",&nc,NULL)); |
57 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose)); |
58 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test BVGetSplit (length %" PetscInt_FMT ", l=%" PetscInt_FMT ", k=%" PetscInt_FMT ", nc=%" PetscInt_FMT ").\n",n,l,k,nc)); |
59 | |||
60 | /* Create template vector */ | ||
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.
|
192 | PetscCall(VecCreate(PETSC_COMM_WORLD,&t)); |
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.
|
192 | PetscCall(VecSetSizes(t,PETSC_DECIDE,n)); |
63 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(VecSetFromOptions(t)); |
64 | |||
65 | /* Create BV object X */ | ||
66 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(BVCreate(PETSC_COMM_WORLD,&X)); |
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.
|
192 | PetscCall(PetscObjectSetName((PetscObject)X,"X")); |
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.
|
192 | PetscCall(BVSetSizesFromVec(X,t,k)); |
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.
|
192 | PetscCall(BVSetFromOptions(X)); |
70 | |||
71 | /* Generate constraints and attach them to X */ | ||
72 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
192 | if (nc>0) { |
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.
|
96 | PetscCall(VecDuplicateVecs(t,nc,&C)); |
74 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
288 | for (j=0;j<nc;j++) { |
75 |
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.
|
480 | for (i=0;i<=j;i++) PetscCall(VecSetValue(C[j],nc-i+1,1.0,INSERT_VALUES)); |
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.
|
192 | PetscCall(VecAssemblyBegin(C[j])); |
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.
|
192 | PetscCall(VecAssemblyEnd(C[j])); |
78 | } | ||
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.
|
96 | PetscCall(BVInsertConstraints(X,&nc,C)); |
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.
|
96 | PetscCall(VecDestroyVecs(nc,&C)); |
81 | } | ||
82 | |||
83 | /* Set up viewer */ | ||
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.
|
192 | PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view)); |
85 |
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.
|
192 | if (verbose) PetscCall(PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB)); |
86 | |||
87 | /* Fill X entries */ | ||
88 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1152 | for (j=0;j<k;j++) { |
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.
|
960 | PetscCall(BVGetColumn(X,j,&v)); |
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.
|
960 | PetscCall(VecSet(v,0.0)); |
91 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
4800 | for (i=0;i<4;i++) { |
92 |
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.
|
3840 | if (i+j<n) PetscCall(VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES)); |
93 | } | ||
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.
|
960 | PetscCall(VecAssemblyBegin(v)); |
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.
|
960 | PetscCall(VecAssemblyEnd(v)); |
96 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
960 | PetscCall(BVRestoreColumn(X,j,&v)); |
97 | } | ||
98 |
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.
|
192 | if (verbose) PetscCall(BVView(X,view)); |
99 | |||
100 | /* Get split BVs */ | ||
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.
|
192 | PetscCall(BVSetActiveColumns(X,l,k)); |
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.
|
192 | PetscCall(BVGetSplit(X,&L,&R)); |
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.
|
192 | PetscCall(PetscObjectSetName((PetscObject)L,"L")); |
104 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(PetscObjectSetName((PetscObject)R,"R")); |
105 | |||
106 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
192 | if (verbose) { |
107 | ✗ | PetscCall(BVView(L,view)); | |
108 | ✗ | PetscCall(BVView(R,view)); | |
109 | } | ||
110 | |||
111 | /* Modify first column of R */ | ||
112 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(BVGetColumn(R,0,&v)); |
113 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(VecSet(v,-1.0)); |
114 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(BVRestoreColumn(R,0,&v)); |
115 | |||
116 | /* Finished using the split BVs */ | ||
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.
|
192 | PetscCall(BVRestoreSplit(X,&L,&R)); |
118 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(PrintFirstRow(X)); |
119 |
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.
|
192 | if (verbose) PetscCall(BVView(X,view)); |
120 | |||
121 | /* Get the left split BV only */ | ||
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.
|
192 | PetscCall(BVGetSplit(X,&L,NULL)); |
123 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
768 | for (j=0;j<l;j++) { |
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.
|
576 | PetscCall(BVOrthogonalizeColumn(L,j,NULL,&norm,NULL)); |
125 | 576 | alpha = 1.0/norm; | |
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.
|
576 | PetscCall(BVScaleColumn(L,j,alpha)); |
127 | } | ||
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.
|
192 | PetscCall(BVRestoreSplit(X,&L,NULL)); |
129 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(PrintFirstRow(X)); |
130 |
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.
|
192 | if (verbose) PetscCall(BVView(X,view)); |
131 | |||
132 | /* Now get the right split BV after changing the number of leading columns */ | ||
133 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(BVSetActiveColumns(X,l-1,k)); |
134 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(BVGetSplit(X,NULL,&R)); |
135 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(BVGetColumn(R,0,&v)); |
136 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(BVInsertVec(X,0,v)); |
137 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(BVRestoreColumn(R,0,&v)); |
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.
|
192 | PetscCall(BVRestoreSplit(X,NULL,&R)); |
139 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
192 | PetscCall(PrintFirstRow(X)); |
140 |
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.
|
192 | if (verbose) PetscCall(BVView(X,view)); |
141 | |||
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.
|
192 | PetscCall(BVDestroy(&X)); |
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.
|
192 | PetscCall(VecDestroy(&t)); |
144 |
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.
|
192 | PetscCall(SlepcFinalize()); |
145 | return 0; | ||
146 | } | ||
147 | |||
148 | /*TEST | ||
149 | |||
150 | testset: | ||
151 | nsize: 2 | ||
152 | output_file: output/test15_1.out | ||
153 | test: | ||
154 | suffix: 1 | ||
155 | args: -bv_type {{vecs contiguous svec mat}shared output} | ||
156 | test: | ||
157 | suffix: 1_cuda | ||
158 | args: -bv_type {{svec mat}} -vec_type cuda | ||
159 | requires: cuda | ||
160 | test: | ||
161 | suffix: 1_hip | ||
162 | args: -bv_type {{svec mat}} -vec_type hip | ||
163 | requires: hip | ||
164 | |||
165 | testset: | ||
166 | nsize: 2 | ||
167 | output_file: output/test15_2.out | ||
168 | test: | ||
169 | suffix: 2 | ||
170 | args: -nc 2 -bv_type {{vecs contiguous svec mat}shared output} | ||
171 | test: | ||
172 | suffix: 2_cuda | ||
173 | args: -nc 2 -bv_type {{svec mat}} -vec_type cuda | ||
174 | requires: cuda | ||
175 | test: | ||
176 | suffix: 2_hip | ||
177 | args: -nc 2 -bv_type {{svec mat}} -vec_type hip | ||
178 | requires: hip | ||
179 | |||
180 | TEST*/ | ||
181 |