| 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 various ST interface functions.\n\n"; | ||
| 12 | |||
| 13 | #include <slepceps.h> | ||
| 14 | |||
| 15 | 45 | int main (int argc,char **argv) | |
| 16 | { | ||
| 17 | 45 | ST st; | |
| 18 | 45 | KSP ksp; | |
| 19 | 45 | PC pc; | |
| 20 | 45 | Mat A,mat[1],Op; | |
| 21 | 45 | Vec v,w; | |
| 22 | 45 | PetscInt N,n=4,i,j,II,Istart,Iend; | |
| 23 | 45 | PetscScalar d,val1,val2; | |
| 24 | 45 | PetscBool test_compl=PETSC_FALSE; | |
| 25 | |||
| 26 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
45 | 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.
|
45 | 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.
|
45 | PetscCall(PetscOptionsGetBool(NULL,NULL,"-complex",&test_compl,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.
|
45 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
| 30 | 45 | N = n*n; | |
| 31 | |||
| 32 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 33 | Create non-symmetric matrix | ||
| 34 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 35 | |||
| 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.
|
45 | PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); |
| 37 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
45 | PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N)); |
| 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.
|
45 | PetscCall(MatSetFromOptions(A)); |
| 39 | |||
| 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.
|
45 | PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); |
| 41 | #if defined(PETSC_USE_COMPLEX) | ||
| 42 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
30 | val1 = test_compl? PetscCMPLX(1.0,-0.2): 1.0; |
| 43 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
30 | val2 = test_compl? PetscCMPLX(-1.0,0.4): -1.0; |
| 44 | #else | ||
| 45 | 15 | val1 = 1.0; | |
| 46 | 15 | val2 = -1.0; | |
| 47 | #endif | ||
| 48 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
765 | for (II=Istart;II<Iend;II++) { |
| 49 | 720 | i = II/n; j = II-i*n; | |
| 50 | 720 | d = 0.0; | |
| 51 |
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.
|
720 | if (i>0) { PetscCall(MatSetValue(A,II,II-n,val1,INSERT_VALUES)); d=d+1.0; } |
| 52 |
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.
|
720 | if (i<n-1) { PetscCall(MatSetValue(A,II,II+n,val2,INSERT_VALUES)); d=d+1.0; } |
| 53 |
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.
|
720 | if (j>0) { PetscCall(MatSetValue(A,II,II-1,val1,INSERT_VALUES)); d=d+1.0; } |
| 54 |
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.
|
720 | if (j<n-1) { PetscCall(MatSetValue(A,II,II+1,val2,INSERT_VALUES)); d=d+1.0; } |
| 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.
|
720 | PetscCall(MatSetValue(A,II,II,d,INSERT_VALUES)); |
| 56 | } | ||
| 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.
|
45 | PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); |
| 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.
|
45 | PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); |
| 59 | |||
| 60 | #if defined(PETSC_USE_COMPLEX) | ||
| 61 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
30 | val1 = test_compl? PetscCMPLX(-0.5,0.01): -0.5; |
| 62 | #else | ||
| 63 | 15 | val1 = -0.5; | |
| 64 | #endif | ||
| 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.
|
45 | PetscCall(MatCreateVecs(A,&v,&w)); |
| 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.
|
45 | PetscCall(VecSetValue(v,0,val1,INSERT_VALUES)); |
| 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.
|
45 | PetscCall(VecSetValue(v,1,1.5,INSERT_VALUES)); |
| 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.
|
45 | PetscCall(VecSetValue(v,2,2,INSERT_VALUES)); |
| 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.
|
45 | PetscCall(VecAssemblyBegin(v)); |
| 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.
|
45 | PetscCall(VecAssemblyEnd(v)); |
| 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.
|
45 | PetscCall(VecView(v,NULL)); |
| 72 | |||
| 73 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 74 | Create the spectral transformation object | ||
| 75 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 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.
|
45 | PetscCall(STCreate(PETSC_COMM_WORLD,&st)); |
| 77 | 45 | mat[0] = A; | |
| 78 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
45 | PetscCall(STSetMatrices(st,1,mat)); |
| 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.
|
45 | PetscCall(STSetType(st,STCAYLEY)); |
| 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.
|
45 | PetscCall(STSetShift(st,2.0)); |
| 81 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
45 | PetscCall(STCayleySetAntishift(st,1.0)); |
| 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.
|
45 | PetscCall(STSetTransform(st,PETSC_TRUE)); |
| 83 | |||
| 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.
|
45 | PetscCall(KSPCreate(PETSC_COMM_WORLD,&ksp)); |
| 85 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
45 | PetscCall(KSPSetType(ksp,KSPPREONLY)); |
| 86 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
45 | PetscCall(KSPGetPC(ksp,&pc)); |
| 87 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
45 | PetscCall(PCSetType(pc,PCLU)); |
| 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.
|
45 | PetscCall(KSPSetTolerances(ksp,100*PETSC_MACHINE_EPSILON,PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT)); |
| 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.
|
45 | PetscCall(KSPSetFromOptions(ksp)); |
| 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.
|
45 | PetscCall(STSetKSP(st,ksp)); |
| 91 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
45 | PetscCall(KSPDestroy(&ksp)); |
| 92 | |||
| 93 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
45 | PetscCall(STSetFromOptions(st)); |
| 94 | |||
| 95 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 96 | Apply the operator | ||
| 97 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 98 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
45 | PetscCall(STApply(st,v,w)); |
| 99 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
45 | PetscCall(VecView(w,NULL)); |
| 100 | |||
| 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.
|
45 | PetscCall(STApplyTranspose(st,v,w)); |
| 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.
|
45 | PetscCall(VecView(w,NULL)); |
| 103 | |||
| 104 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 10 times.
|
45 | if (test_compl) { |
| 105 |
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.
|
15 | PetscCall(STApplyHermitianTranspose(st,v,w)); |
| 106 |
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.
|
15 | PetscCall(VecView(w,NULL)); |
| 107 | } | ||
| 108 | |||
| 109 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
45 | PetscCall(STMatMult(st,1,v,w)); |
| 110 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
45 | PetscCall(VecView(w,NULL)); |
| 111 | |||
| 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.
|
45 | PetscCall(STMatMultTranspose(st,1,v,w)); |
| 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.
|
45 | PetscCall(VecView(w,NULL)); |
| 114 | |||
| 115 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 10 times.
|
45 | if (test_compl) { |
| 116 |
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.
|
15 | PetscCall(STMatMultHermitianTranspose(st,1,v,w)); |
| 117 |
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.
|
15 | PetscCall(VecView(w,NULL)); |
| 118 | } | ||
| 119 | |||
| 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.
|
45 | PetscCall(STMatSolve(st,v,w)); |
| 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.
|
45 | PetscCall(VecView(w,NULL)); |
| 122 | |||
| 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.
|
45 | PetscCall(STMatSolveTranspose(st,v,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.
|
45 | PetscCall(VecView(w,NULL)); |
| 125 | |||
| 126 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 10 times.
|
45 | if (test_compl) { |
| 127 |
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.
|
15 | PetscCall(STMatSolveHermitianTranspose(st,v,w)); |
| 128 |
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.
|
15 | PetscCall(VecView(w,NULL)); |
| 129 | } | ||
| 130 | |||
| 131 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 132 | Get the operator matrix | ||
| 133 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 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.
|
45 | PetscCall(STGetOperator(st,&Op)); |
| 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.
|
45 | PetscCall(MatMult(Op,v,w)); |
| 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.
|
45 | PetscCall(VecView(w,NULL)); |
| 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.
|
45 | PetscCall(MatMultTranspose(Op,v,w)); |
| 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.
|
45 | PetscCall(VecView(w,NULL)); |
| 139 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 10 times.
|
45 | if (test_compl) { |
| 140 |
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.
|
15 | PetscCall(MatMultHermitianTranspose(Op,v,w)); |
| 141 |
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.
|
15 | PetscCall(VecView(w,NULL)); |
| 142 | } | ||
| 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.
|
45 | PetscCall(STRestoreOperator(st,&Op)); |
| 144 | |||
| 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.
|
45 | PetscCall(STDestroy(&st)); |
| 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.
|
45 | PetscCall(MatDestroy(&A)); |
| 147 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
45 | PetscCall(VecDestroy(&v)); |
| 148 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
45 | PetscCall(VecDestroy(&w)); |
| 149 |
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.
|
45 | PetscCall(SlepcFinalize()); |
| 150 | return 0; | ||
| 151 | } | ||
| 152 | |||
| 153 | /*TEST | ||
| 154 | |||
| 155 | testset: | ||
| 156 | output_file: output/test5_1.out | ||
| 157 | requires: !single | ||
| 158 | test: | ||
| 159 | suffix: 1 | ||
| 160 | args: -st_matmode {{copy inplace}} | ||
| 161 | test: | ||
| 162 | suffix: 1_shell | ||
| 163 | args: -st_matmode shell -ksp_type bcgs -pc_type jacobi | ||
| 164 | |||
| 165 | testset: | ||
| 166 | output_file: output/test5_2.out | ||
| 167 | requires: complex !single | ||
| 168 | args: -complex | ||
| 169 | test: | ||
| 170 | suffix: 2 | ||
| 171 | args: -st_matmode {{copy inplace}} | ||
| 172 | test: | ||
| 173 | suffix: 2_shell | ||
| 174 | args: -st_matmode shell -ksp_type bcgs -pc_type jacobi | ||
| 175 | |||
| 176 | TEST*/ | ||
| 177 |