| 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 | Phi functions | ||
| 12 | phi_0(x) = exp(x) | ||
| 13 | phi_1(x) = (exp(x)-1)/x | ||
| 14 | phi_k(x) = (phi_{k-1}(x)-1/(k-1)!)/x | ||
| 15 | */ | ||
| 16 | |||
| 17 | #include <slepc/private/fnimpl.h> /*I "slepcfn.h" I*/ | ||
| 18 | |||
| 19 | typedef struct { | ||
| 20 | PetscInt k; /* index of the phi-function, defaults to k=1 */ | ||
| 21 | Mat H; /* auxiliary matrix of order m+k */ | ||
| 22 | Mat F; /* auxiliary matrix to store exp(H) */ | ||
| 23 | } FN_PHI; | ||
| 24 | |||
| 25 | #define MAX_INDEX 10 | ||
| 26 | |||
| 27 | static const PetscReal rfactorial[MAX_INDEX+2] = { 1, 1, 0.5, 1.0/6, 1.0/24, 1.0/120, 1.0/720, 1.0/5040, 1.0/40320, 1.0/362880, 1.0/3628800, 1.0/39916800 }; | ||
| 28 | |||
| 29 | 270 | static PetscErrorCode FNEvaluateFunction_Phi(FN fn,PetscScalar x,PetscScalar *y) | |
| 30 | { | ||
| 31 | 270 | FN_PHI *ctx = (FN_PHI*)fn->data; | |
| 32 | 270 | PetscInt i; | |
| 33 | 270 | PetscScalar phi[MAX_INDEX+1]; | |
| 34 | |||
| 35 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
270 | PetscFunctionBegin; |
| 36 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
270 | if (x==0.0) *y = rfactorial[ctx->k]; |
| 37 | else { | ||
| 38 | 270 | phi[0] = PetscExpScalar(x); | |
| 39 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
630 | for (i=1;i<=ctx->k;i++) phi[i] = (phi[i-1]-rfactorial[i-1])/x; |
| 40 | 270 | *y = phi[ctx->k]; | |
| 41 | } | ||
| 42 |
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.
|
270 | PetscFunctionReturn(PETSC_SUCCESS); |
| 43 | } | ||
| 44 | |||
| 45 | 30 | static PetscErrorCode FNEvaluateDerivative_Phi(FN fn,PetscScalar x,PetscScalar *y) | |
| 46 | { | ||
| 47 | 30 | FN_PHI *ctx = (FN_PHI*)fn->data; | |
| 48 | 30 | PetscInt i; | |
| 49 | 30 | PetscScalar phi[MAX_INDEX+2]; | |
| 50 | |||
| 51 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
30 | PetscFunctionBegin; |
| 52 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
30 | if (x==0.0) *y = rfactorial[ctx->k+1]; |
| 53 | else { | ||
| 54 | 30 | phi[0] = PetscExpScalar(x); | |
| 55 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
100 | for (i=1;i<=ctx->k+1;i++) phi[i] = (phi[i-1]-rfactorial[i-1])/x; |
| 56 | 30 | *y = phi[ctx->k] - phi[ctx->k+1]*(PetscReal)ctx->k; | |
| 57 | } | ||
| 58 |
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.
|
30 | PetscFunctionReturn(PETSC_SUCCESS); |
| 59 | } | ||
| 60 | |||
| 61 | 230 | static PetscErrorCode FNEvaluateFunctionMatVec_Phi(FN fn,Mat A,Vec v) | |
| 62 | { | ||
| 63 | 230 | FN_PHI *ctx = (FN_PHI*)fn->data; | |
| 64 | 230 | PetscInt i,j,m,n,nh; | |
| 65 | 230 | PetscScalar *Ha,*va,sfactor=1.0; | |
| 66 | 230 | const PetscScalar *Aa,*Fa; | |
| 67 | |||
| 68 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
230 | PetscFunctionBegin; |
| 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.
|
230 | PetscCall(MatGetSize(A,&m,NULL)); |
| 70 | 230 | n = m+ctx->k; | |
| 71 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
230 | if (ctx->H) { |
| 72 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
190 | PetscCall(MatGetSize(ctx->H,&nh,NULL)); |
| 73 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
190 | if (n!=nh) { |
| 74 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
190 | PetscCall(MatDestroy(&ctx->H)); |
| 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.
|
190 | PetscCall(MatDestroy(&ctx->F)); |
| 76 | } | ||
| 77 | } | ||
| 78 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
230 | if (!ctx->H) { |
| 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.
|
230 | PetscCall(MatCreateDense(PETSC_COMM_SELF,n,n,n,n,NULL,&ctx->H)); |
| 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.
|
230 | PetscCall(MatCreateDense(PETSC_COMM_SELF,n,n,n,n,NULL,&ctx->F)); |
| 81 | } | ||
| 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.
|
230 | PetscCall(MatDenseGetArray(ctx->H,&Ha)); |
| 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.
|
230 | PetscCall(MatDenseGetArrayRead(A,&Aa)); |
| 84 |
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.
|
15470 | for (j=0;j<m;j++) PetscCall(PetscArraycpy(Ha+j*n,Aa+j*m,m)); |
| 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.
|
230 | PetscCall(MatDenseRestoreArrayRead(A,&Aa)); |
| 86 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
230 | if (ctx->k) { |
| 87 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
30700 | for (j=0;j<m;j++) for (i=m;i<n;i++) Ha[i+j*n] = 0.0; |
| 88 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
16080 | for (j=m;j<n;j++) for (i=0;i<n;i++) Ha[i+j*n] = 0.0; |
| 89 | 220 | Ha[0+m*n] = fn->alpha; | |
| 90 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
240 | for (j=m+1;j<n;j++) Ha[j-1+j*n] = fn->alpha; |
| 91 | } | ||
| 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.
|
230 | PetscCall(MatDenseRestoreArray(ctx->H,&Ha)); |
| 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.
|
230 | PetscCall(FNEvaluateFunctionMat_Exp_Higham(fn,ctx->H,ctx->F)); |
| 95 | |||
| 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.
|
230 | PetscCall(MatDenseGetArrayRead(ctx->F,&Fa)); |
| 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.
|
230 | PetscCall(VecGetArray(v,&va)); |
| 98 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
230 | if (ctx->k) { |
| 99 | 220 | sfactor = PetscPowScalarInt(fn->alpha,-ctx->k); | |
| 100 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
15600 | for (i=0;i<m;i++) va[i] = sfactor*Fa[i+(n-1)*n]; |
| 101 | } else { | ||
| 102 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
90 | for (i=0;i<m;i++) va[i] = sfactor*Fa[i+0*n]; |
| 103 | } | ||
| 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.
|
230 | PetscCall(VecRestoreArray(v,&va)); |
| 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.
|
230 | PetscCall(MatDenseRestoreArrayRead(ctx->F,&Fa)); |
| 106 |
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.
|
46 | PetscFunctionReturn(PETSC_SUCCESS); |
| 107 | } | ||
| 108 | |||
| 109 | 30 | static PetscErrorCode FNPhiSetIndex_Phi(FN fn,PetscInt k) | |
| 110 | { | ||
| 111 | 30 | FN_PHI *ctx = (FN_PHI*)fn->data; | |
| 112 | |||
| 113 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
30 | PetscFunctionBegin; |
| 114 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
30 | PetscCheck(k>=0,PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Index cannot be negative"); |
| 115 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
30 | PetscCheck(k<=MAX_INDEX,PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Phi functions only implemented for k<=%d",MAX_INDEX); |
| 116 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
30 | if (k!=ctx->k) { |
| 117 | 20 | ctx->k = k; | |
| 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.
|
20 | PetscCall(MatDestroy(&ctx->H)); |
| 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.
|
20 | PetscCall(MatDestroy(&ctx->F)); |
| 120 | } | ||
| 121 |
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.
|
6 | PetscFunctionReturn(PETSC_SUCCESS); |
| 122 | } | ||
| 123 | |||
| 124 | /*@ | ||
| 125 | FNPhiSetIndex - Sets the index of the $\varphi$-function. | ||
| 126 | |||
| 127 | Logically Collective | ||
| 128 | |||
| 129 | Input Parameters: | ||
| 130 | + fn - the math function context | ||
| 131 | - k - the index | ||
| 132 | |||
| 133 | Notes: | ||
| 134 | The $\varphi$-functions are defined as | ||
| 135 | $$\begin{align} | ||
| 136 | \varphi_0(x)&=e^x,\\ | ||
| 137 | \varphi_1(x)&=\frac{e^x-1}{x},\\ | ||
| 138 | \varphi_k(x)&=\frac{\varphi_{k-1}(x)-1/(k-1)!}{x}. | ||
| 139 | \end{align}$$ | ||
| 140 | The default index is $k=1$. | ||
| 141 | |||
| 142 | Level: intermediate | ||
| 143 | |||
| 144 | .seealso: [](sec:fn), `FNPHI`, `FNPhiGetIndex()` | ||
| 145 | @*/ | ||
| 146 | 30 | PetscErrorCode FNPhiSetIndex(FN fn,PetscInt k) | |
| 147 | { | ||
| 148 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
30 | PetscFunctionBegin; |
| 149 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
30 | PetscValidHeaderSpecific(fn,FN_CLASSID,1); |
| 150 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ 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 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
30 | PetscValidLogicalCollectiveInt(fn,k,2); |
| 151 |
8/14✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 10 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.
|
30 | PetscTryMethod(fn,"FNPhiSetIndex_C",(FN,PetscInt),(fn,k)); |
| 152 |
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.
|
30 | PetscFunctionReturn(PETSC_SUCCESS); |
| 153 | } | ||
| 154 | |||
| 155 | 10 | static PetscErrorCode FNPhiGetIndex_Phi(FN fn,PetscInt *k) | |
| 156 | { | ||
| 157 | 10 | FN_PHI *ctx = (FN_PHI*)fn->data; | |
| 158 | |||
| 159 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
| 160 | 10 | *k = ctx->k; | |
| 161 |
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.
|
10 | PetscFunctionReturn(PETSC_SUCCESS); |
| 162 | } | ||
| 163 | |||
| 164 | /*@ | ||
| 165 | FNPhiGetIndex - Gets the index of the $\varphi$-function. | ||
| 166 | |||
| 167 | Not Collective | ||
| 168 | |||
| 169 | Input Parameter: | ||
| 170 | . fn - the math function context | ||
| 171 | |||
| 172 | Output Parameter: | ||
| 173 | . k - the index | ||
| 174 | |||
| 175 | Level: intermediate | ||
| 176 | |||
| 177 | .seealso: [](sec:fn), `FNPHI`, `FNPhiSetIndex()` | ||
| 178 | @*/ | ||
| 179 | 10 | PetscErrorCode FNPhiGetIndex(FN fn,PetscInt *k) | |
| 180 | { | ||
| 181 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
| 182 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
10 | PetscValidHeaderSpecific(fn,FN_CLASSID,1); |
| 183 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
10 | PetscAssertPointer(k,2); |
| 184 |
9/16✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 10 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
|
10 | PetscUseMethod(fn,"FNPhiGetIndex_C",(FN,PetscInt*),(fn,k)); |
| 185 |
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.
|
10 | PetscFunctionReturn(PETSC_SUCCESS); |
| 186 | } | ||
| 187 | |||
| 188 | 30 | static PetscErrorCode FNView_Phi(FN fn,PetscViewer viewer) | |
| 189 | { | ||
| 190 | 30 | FN_PHI *ctx = (FN_PHI*)fn->data; | |
| 191 | 30 | PetscBool isascii; | |
| 192 | 30 | char str[50],strx[50]; | |
| 193 | |||
| 194 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
30 | PetscFunctionBegin; |
| 195 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); |
| 196 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
30 | if (isascii) { |
| 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.
|
30 | PetscCall(PetscViewerASCIIPrintf(viewer," phi_%" PetscInt_FMT ": ",ctx->k)); |
| 198 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); |
| 199 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
30 | if (fn->beta!=(PetscScalar)1.0) { |
| 200 |
4/6✓ Branch 0 taken 2 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(SlepcSNPrintfScalar(str,sizeof(str),fn->beta,PETSC_TRUE)); |
| 201 |
4/6✓ Branch 0 taken 2 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(PetscViewerASCIIPrintf(viewer,"%s*",str)); |
| 202 | } | ||
| 203 |
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.
|
30 | if (fn->alpha==(PetscScalar)1.0) PetscCall(PetscSNPrintf(strx,sizeof(strx),"x")); |
| 204 | else { | ||
| 205 |
4/6✓ Branch 0 taken 2 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(SlepcSNPrintfScalar(str,sizeof(str),fn->alpha,PETSC_TRUE)); |
| 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.
|
10 | PetscCall(PetscSNPrintf(strx,sizeof(strx),"(%s*x)",str)); |
| 207 | } | ||
| 208 |
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.
|
30 | if (!ctx->k) PetscCall(PetscViewerASCIIPrintf(viewer,"exp(%s)\n",strx)); |
| 209 |
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.
|
20 | else if (ctx->k==1) PetscCall(PetscViewerASCIIPrintf(viewer,"(exp(%s)-1)/%s\n",strx,strx)); |
| 210 |
4/6✓ Branch 0 taken 2 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 | else PetscCall(PetscViewerASCIIPrintf(viewer,"(phi_%" PetscInt_FMT "(%s)-1/%" PetscInt_FMT "!)/%s\n",ctx->k-1,strx,ctx->k-1,strx)); |
| 211 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); |
| 212 | } | ||
| 213 |
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.
|
6 | PetscFunctionReturn(PETSC_SUCCESS); |
| 214 | } | ||
| 215 | |||
| 216 | 20 | static PetscErrorCode FNSetFromOptions_Phi(FN fn,PetscOptionItems PetscOptionsObject) | |
| 217 | { | ||
| 218 | 20 | FN_PHI *ctx = (FN_PHI*)fn->data; | |
| 219 | 20 | PetscInt k; | |
| 220 | 20 | PetscBool flag; | |
| 221 | |||
| 222 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
20 | PetscFunctionBegin; |
| 223 |
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.
|
20 | PetscOptionsHeadBegin(PetscOptionsObject,"FN Phi Options"); |
| 224 | |||
| 225 |
4/6✓ Branch 0 taken 2 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(PetscOptionsInt("-fn_phi_index","Index of the phi-function","FNPhiSetIndex",ctx->k,&k,&flag)); |
| 226 |
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.
|
20 | if (flag) PetscCall(FNPhiSetIndex(fn,k)); |
| 227 | |||
| 228 |
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.
|
20 | PetscOptionsHeadEnd(); |
| 229 |
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); |
| 230 | } | ||
| 231 | |||
| 232 | 10 | static PetscErrorCode FNDuplicate_Phi(FN fn,MPI_Comm comm,FN *newfn) | |
| 233 | { | ||
| 234 | 10 | FN_PHI *ctx = (FN_PHI*)fn->data,*ctx2 = (FN_PHI*)(*newfn)->data; | |
| 235 | |||
| 236 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
| 237 | 10 | ctx2->k = ctx->k; | |
| 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.
|
10 | PetscFunctionReturn(PETSC_SUCCESS); |
| 239 | } | ||
| 240 | |||
| 241 | 50 | static PetscErrorCode FNDestroy_Phi(FN fn) | |
| 242 | { | ||
| 243 | 50 | FN_PHI *ctx = (FN_PHI*)fn->data; | |
| 244 | |||
| 245 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
50 | PetscFunctionBegin; |
| 246 |
4/6✓ Branch 0 taken 2 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(MatDestroy(&ctx->H)); |
| 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.
|
50 | PetscCall(MatDestroy(&ctx->F)); |
| 248 |
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.
|
50 | PetscCall(PetscFree(fn->data)); |
| 249 |
4/6✓ Branch 0 taken 2 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(PetscObjectComposeFunction((PetscObject)fn,"FNPhiSetIndex_C",NULL)); |
| 250 |
4/6✓ Branch 0 taken 2 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(PetscObjectComposeFunction((PetscObject)fn,"FNPhiGetIndex_C",NULL)); |
| 251 |
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.
|
10 | PetscFunctionReturn(PETSC_SUCCESS); |
| 252 | } | ||
| 253 | |||
| 254 | /*MC | ||
| 255 | FNPHI - FNPHI = "phi" - One of the $\varphi$-functions $\varphi_0(x)$, $\varphi_1(x)$, ..., | ||
| 256 | where the index of the $\varphi$-function is specified with `FNPhiSetIndex()`. | ||
| 257 | |||
| 258 | Level: beginner | ||
| 259 | |||
| 260 | Notes: | ||
| 261 | The $\varphi$-functions are defined as | ||
| 262 | $$\begin{align} | ||
| 263 | \varphi_0(x)&=e^x,\\ | ||
| 264 | \varphi_1(x)&=\frac{e^x-1}{x},\\ | ||
| 265 | \varphi_k(x)&=\frac{\varphi_{k-1}(x)-1/(k-1)!}{x}. | ||
| 266 | \end{align}$$ | ||
| 267 | The default index is $k=1$. | ||
| 268 | |||
| 269 | .seealso: [](sec:fn), `FN`, `FNType`, `FNSetType()`, `FNPhiSetIndex()` | ||
| 270 | M*/ | ||
| 271 | |||
| 272 | 50 | SLEPC_EXTERN PetscErrorCode FNCreate_Phi(FN fn) | |
| 273 | { | ||
| 274 | 50 | FN_PHI *ctx; | |
| 275 | |||
| 276 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
50 | PetscFunctionBegin; |
| 277 |
4/6✓ Branch 0 taken 2 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(PetscNew(&ctx)); |
| 278 | 50 | fn->data = (void*)ctx; | |
| 279 | 50 | ctx->k = 1; | |
| 280 | |||
| 281 | 50 | fn->ops->evaluatefunction = FNEvaluateFunction_Phi; | |
| 282 | 50 | fn->ops->evaluatederivative = FNEvaluateDerivative_Phi; | |
| 283 | 50 | fn->ops->evaluatefunctionmatvec[0] = FNEvaluateFunctionMatVec_Phi; | |
| 284 | 50 | fn->ops->setfromoptions = FNSetFromOptions_Phi; | |
| 285 | 50 | fn->ops->view = FNView_Phi; | |
| 286 | 50 | fn->ops->duplicate = FNDuplicate_Phi; | |
| 287 | 50 | fn->ops->destroy = FNDestroy_Phi; | |
| 288 |
4/6✓ Branch 0 taken 2 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(PetscObjectComposeFunction((PetscObject)fn,"FNPhiSetIndex_C",FNPhiSetIndex_Phi)); |
| 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(PetscObjectComposeFunction((PetscObject)fn,"FNPhiGetIndex_C",FNPhiGetIndex_Phi)); |
| 290 |
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.
|
10 | PetscFunctionReturn(PETSC_SUCCESS); |
| 291 | } | ||
| 292 |