| 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 matrix exponential.\n\n"; | ||
| 12 | |||
| 13 | #include <slepcfn.h> | ||
| 14 | |||
| 15 | /* | ||
| 16 | Compute matrix exponential B = expm(A) | ||
| 17 | */ | ||
| 18 | 432 | PetscErrorCode TestMatExp(FN fn,Mat A,PetscViewer viewer,PetscBool verbose,PetscBool inplace,PetscBool checkerror) | |
| 19 | { | ||
| 20 | 432 | PetscScalar tau,eta; | |
| 21 | 432 | PetscBool set,flg; | |
| 22 | 432 | PetscInt n; | |
| 23 | 432 | Mat F,R,Finv,Acopy; | |
| 24 | 432 | Vec v,f0; | |
| 25 | 432 | FN finv; | |
| 26 | 432 | PetscReal nrm,nrmf; | |
| 27 | |||
| 28 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
432 | PetscFunctionBeginUser; |
| 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.
|
432 | PetscCall(MatGetSize(A,&n,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.
|
432 | PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&F)); |
| 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.
|
432 | PetscCall(PetscObjectSetName((PetscObject)F,"F")); |
| 32 | /* compute matrix exponential */ | ||
| 33 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
432 | if (inplace) { |
| 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.
|
96 | PetscCall(MatCopy(A,F,SAME_NONZERO_PATTERN)); |
| 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.
|
96 | PetscCall(MatIsHermitianKnown(A,&set,&flg)); |
| 36 |
7/10✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
96 | if (set && flg) PetscCall(MatSetOption(F,MAT_HERMITIAN,PETSC_TRUE)); |
| 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.
|
96 | PetscCall(FNEvaluateFunctionMat(fn,F,NULL)); |
| 38 | } else { | ||
| 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.
|
336 | PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&Acopy)); |
| 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.
|
336 | PetscCall(FNEvaluateFunctionMat(fn,A,F)); |
| 41 | /* check that A has not been modified */ | ||
| 42 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
336 | PetscCall(MatAXPY(Acopy,-1.0,A,SAME_NONZERO_PATTERN)); |
| 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.
|
336 | PetscCall(MatNorm(Acopy,NORM_1,&nrm)); |
| 44 |
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.
|
336 | if (nrm>100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: the input matrix has changed by %g\n",(double)nrm)); |
| 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.
|
336 | PetscCall(MatDestroy(&Acopy)); |
| 46 | } | ||
| 47 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
432 | if (verbose) { |
| 48 | ✗ | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n")); | |
| 49 | ✗ | PetscCall(MatView(A,viewer)); | |
| 50 | ✗ | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed expm(A) - - - - - - -\n")); | |
| 51 | ✗ | PetscCall(MatView(F,viewer)); | |
| 52 | } | ||
| 53 | /* print matrix norm for checking */ | ||
| 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.
|
432 | PetscCall(MatNorm(F,NORM_1,&nrmf)); |
| 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.
|
432 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"The 1-norm of f(A) is %g\n",(double)nrmf)); |
| 56 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
432 | if (checkerror) { |
| 57 | ✗ | PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Finv)); | |
| 58 | ✗ | PetscCall(PetscObjectSetName((PetscObject)Finv,"Finv")); | |
| 59 | ✗ | PetscCall(FNGetScale(fn,&tau,&eta)); | |
| 60 | /* compute inverse exp(-tau*A)/eta */ | ||
| 61 | ✗ | PetscCall(FNCreate(PETSC_COMM_WORLD,&finv)); | |
| 62 | ✗ | PetscCall(FNSetType(finv,FNEXP)); | |
| 63 | ✗ | PetscCall(FNSetFromOptions(finv)); | |
| 64 | ✗ | PetscCall(FNSetScale(finv,-tau,1.0/eta)); | |
| 65 | ✗ | if (inplace) { | |
| 66 | ✗ | PetscCall(MatCopy(A,Finv,SAME_NONZERO_PATTERN)); | |
| 67 | ✗ | PetscCall(MatIsHermitianKnown(A,&set,&flg)); | |
| 68 | ✗ | if (set && flg) PetscCall(MatSetOption(Finv,MAT_HERMITIAN,PETSC_TRUE)); | |
| 69 | ✗ | PetscCall(FNEvaluateFunctionMat(finv,Finv,NULL)); | |
| 70 | ✗ | } else PetscCall(FNEvaluateFunctionMat(finv,A,Finv)); | |
| 71 | ✗ | PetscCall(FNDestroy(&finv)); | |
| 72 | /* check error ||F*Finv-I||_F */ | ||
| 73 | ✗ | PetscCall(MatMatMult(F,Finv,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&R)); | |
| 74 | ✗ | PetscCall(MatShift(R,-1.0)); | |
| 75 | ✗ | PetscCall(MatNorm(R,NORM_FROBENIUS,&nrm)); | |
| 76 | ✗ | if (nrm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"||exp(A)*exp(-A)-I||_F < 100*eps\n")); | |
| 77 | ✗ | else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"||exp(A)*exp(-A)-I||_F = %g\n",(double)nrm)); | |
| 78 | ✗ | PetscCall(MatDestroy(&R)); | |
| 79 | ✗ | PetscCall(MatDestroy(&Finv)); | |
| 80 | } | ||
| 81 | /* check FNEvaluateFunctionMatVec() */ | ||
| 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.
|
432 | PetscCall(MatCreateVecs(A,&v,&f0)); |
| 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.
|
432 | PetscCall(MatGetColumnVector(F,f0,0)); |
| 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.
|
432 | PetscCall(FNEvaluateFunctionMatVec(fn,A,v)); |
| 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.
|
432 | PetscCall(VecAXPY(v,-1.0,f0)); |
| 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.
|
432 | PetscCall(VecNorm(v,NORM_2,&nrm)); |
| 87 |
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.
|
432 | if (nrm/nrmf>100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of f(A)*e_1-v is %g\n",(double)nrm)); |
| 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.
|
432 | PetscCall(MatDestroy(&F)); |
| 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.
|
432 | PetscCall(VecDestroy(&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.
|
432 | PetscCall(VecDestroy(&f0)); |
| 91 |
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.
|
80 | PetscFunctionReturn(PETSC_SUCCESS); |
| 92 | } | ||
| 93 | |||
| 94 | 216 | int main(int argc,char **argv) | |
| 95 | { | ||
| 96 | 216 | FN fn; | |
| 97 | 216 | Mat A=NULL; | |
| 98 | 216 | PetscInt i,j,n=10; | |
| 99 | 216 | PetscScalar *As; | |
| 100 | 216 | PetscViewer viewer; | |
| 101 | 216 | PetscBool verbose,inplace,checkerror,matcuda; | |
| 102 | |||
| 103 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
216 | PetscFunctionBeginUser; |
| 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.
|
216 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
| 105 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
216 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
| 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.
|
216 | PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose)); |
| 107 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
216 | PetscCall(PetscOptionsHasName(NULL,NULL,"-inplace",&inplace)); |
| 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.
|
216 | PetscCall(PetscOptionsHasName(NULL,NULL,"-checkerror",&checkerror)); |
| 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.
|
216 | PetscCall(PetscOptionsHasName(NULL,NULL,"-matcuda",&matcuda)); |
| 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.
|
216 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Matrix exponential, n=%" PetscInt_FMT ".\n",n)); |
| 111 | |||
| 112 | /* Create exponential function object */ | ||
| 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.
|
216 | PetscCall(FNCreate(PETSC_COMM_WORLD,&fn)); |
| 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.
|
216 | PetscCall(FNSetType(fn,FNEXP)); |
| 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.
|
216 | PetscCall(FNSetFromOptions(fn)); |
| 116 | |||
| 117 | /* Set up viewer */ | ||
| 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.
|
216 | PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer)); |
| 119 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
216 | PetscCall(FNView(fn,viewer)); |
| 120 |
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.
|
216 | if (verbose) PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB)); |
| 121 | |||
| 122 | /* Create matrices */ | ||
| 123 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 2 times.
|
216 | if (matcuda) { |
| 124 | #if defined(PETSC_HAVE_CUDA) | ||
| 125 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
16 | PetscCall(MatCreateSeqDenseCUDA(PETSC_COMM_SELF,n,n,NULL,&A)); |
| 126 | #endif | ||
| 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.
|
200 | } else PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A)); |
| 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.
|
216 | PetscCall(PetscObjectSetName((PetscObject)A,"A")); |
| 129 | |||
| 130 | /* Fill A with a symmetric Toeplitz matrix */ | ||
| 131 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
216 | PetscCall(MatDenseGetArray(A,&As)); |
| 132 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
6776 | for (i=0;i<n;i++) As[i+i*n]=2.0; |
| 133 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
648 | for (j=1;j<3;j++) { |
| 134 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
12904 | for (i=0;i<n-j;i++) { |
| 135 | 12472 | As[i+(i+j)*n]=1.0; | |
| 136 | 12472 | As[(i+j)+i*n]=1.0; | |
| 137 | } | ||
| 138 | } | ||
| 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.
|
216 | PetscCall(MatDenseRestoreArray(A,&As)); |
| 140 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
216 | PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE)); |
| 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.
|
216 | PetscCall(TestMatExp(fn,A,viewer,verbose,inplace,checkerror)); |
| 142 | |||
| 143 | /* Repeat with non-symmetric A */ | ||
| 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.
|
216 | PetscCall(MatDenseGetArray(A,&As)); |
| 145 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
648 | for (j=1;j<3;j++) { |
| 146 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
12904 | for (i=0;i<n-j;i++) As[(i+j)+i*n]=-1.0; |
| 147 | } | ||
| 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.
|
216 | PetscCall(MatDenseRestoreArray(A,&As)); |
| 149 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
216 | PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE)); |
| 150 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
216 | PetscCall(TestMatExp(fn,A,viewer,verbose,inplace,checkerror)); |
| 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.
|
216 | PetscCall(MatDestroy(&A)); |
| 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.
|
216 | PetscCall(FNDestroy(&fn)); |
| 154 |
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.
|
216 | PetscCall(SlepcFinalize()); |
| 155 | return 0; | ||
| 156 | } | ||
| 157 | |||
| 158 | /*TEST | ||
| 159 | |||
| 160 | testset: | ||
| 161 | filter: grep -v "computing matrix functions" | ||
| 162 | output_file: output/test3_1.out | ||
| 163 | test: | ||
| 164 | suffix: 1 | ||
| 165 | args: -fn_method {{0 1}} | ||
| 166 | test: | ||
| 167 | suffix: 1_subdiagonalpade | ||
| 168 | args: -fn_method {{2 3}} | ||
| 169 | requires: c99_complex !single | ||
| 170 | test: | ||
| 171 | suffix: 1_cuda | ||
| 172 | args: -fn_method 1 -matcuda | ||
| 173 | requires: cuda !magma | ||
| 174 | test: | ||
| 175 | suffix: 1_magma | ||
| 176 | args: -fn_method {{0 1 2 3}} -matcuda | ||
| 177 | requires: cuda magma | ||
| 178 | test: | ||
| 179 | suffix: 2 | ||
| 180 | args: -inplace -fn_method{{0 1}} | ||
| 181 | test: | ||
| 182 | suffix: 2_subdiagonalpade | ||
| 183 | args: -inplace -fn_method{{2 3}} | ||
| 184 | requires: c99_complex !single | ||
| 185 | test: | ||
| 186 | suffix: 2_cuda | ||
| 187 | args: -inplace -fn_method 1 -matcuda | ||
| 188 | requires: cuda !magma | ||
| 189 | test: | ||
| 190 | suffix: 2_magma | ||
| 191 | args: -inplace -fn_method {{0 1 2 3}} -matcuda | ||
| 192 | requires: cuda magma | ||
| 193 | |||
| 194 | testset: | ||
| 195 | args: -fn_scale 0.1 | ||
| 196 | filter: grep -v "computing matrix functions" | ||
| 197 | output_file: output/test3_3.out | ||
| 198 | test: | ||
| 199 | suffix: 3 | ||
| 200 | args: -fn_method {{0 1}} | ||
| 201 | test: | ||
| 202 | suffix: 3_subdiagonalpade | ||
| 203 | args: -fn_method {{2 3}} | ||
| 204 | requires: c99_complex !single | ||
| 205 | |||
| 206 | testset: | ||
| 207 | args: -n 120 -fn_scale 0.6,1.5 | ||
| 208 | filter: grep -v "computing matrix functions" | ||
| 209 | output_file: output/test3_4.out | ||
| 210 | test: | ||
| 211 | suffix: 4 | ||
| 212 | args: -fn_method {{0 1}} | ||
| 213 | requires: !single | ||
| 214 | test: | ||
| 215 | suffix: 4_subdiagonalpade | ||
| 216 | args: -fn_method {{2 3}} | ||
| 217 | requires: c99_complex !single | ||
| 218 | |||
| 219 | test: | ||
| 220 | suffix: 5 | ||
| 221 | args: -fn_scale 30 -fn_method {{2 3}} | ||
| 222 | filter: grep -v "computing matrix functions" | ||
| 223 | requires: c99_complex !single | ||
| 224 | output_file: output/test3_5.out | ||
| 225 | |||
| 226 | test: | ||
| 227 | suffix: 6 | ||
| 228 | args: -fn_scale 1e-9 -fn_method {{2 3}} | ||
| 229 | filter: grep -v "computing matrix functions" | ||
| 230 | requires: c99_complex !single | ||
| 231 | output_file: output/test3_6.out | ||
| 232 | |||
| 233 | TEST*/ | ||
| 234 |