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 | Rational function r(x) = p(x)/q(x), where p(x) and q(x) are polynomials | ||
12 | */ | ||
13 | |||
14 | #include <slepc/private/fnimpl.h> /*I "slepcfn.h" I*/ | ||
15 | |||
16 | typedef struct { | ||
17 | PetscScalar *pcoeff; /* numerator coefficients */ | ||
18 | PetscInt np; /* length of array pcoeff, p(x) has degree np-1 */ | ||
19 | PetscScalar *qcoeff; /* denominator coefficients */ | ||
20 | PetscInt nq; /* length of array qcoeff, q(x) has degree nq-1 */ | ||
21 | } FN_RATIONAL; | ||
22 | |||
23 | 93511 | static PetscErrorCode FNEvaluateFunction_Rational(FN fn,PetscScalar x,PetscScalar *y) | |
24 | { | ||
25 | 93511 | FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data; | |
26 | 93511 | PetscInt i; | |
27 | 93511 | PetscScalar p,q; | |
28 | |||
29 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
93511 | PetscFunctionBegin; |
30 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
93511 | if (!ctx->np) p = 1.0; |
31 | else { | ||
32 | 93495 | p = ctx->pcoeff[0]; | |
33 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
152124 | for (i=1;i<ctx->np;i++) |
34 | 58629 | p = ctx->pcoeff[i]+x*p; | |
35 | } | ||
36 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
93511 | if (!ctx->nq) *y = p; |
37 | else { | ||
38 | 12067 | q = ctx->qcoeff[0]; | |
39 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
25038 | for (i=1;i<ctx->nq;i++) |
40 | 12971 | q = ctx->qcoeff[i]+x*q; | |
41 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
12067 | PetscCheck(q!=0.0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Function not defined in the requested value"); |
42 | 12067 | *y = p/q; | |
43 | } | ||
44 |
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.
|
18536 | PetscFunctionReturn(PETSC_SUCCESS); |
45 | } | ||
46 | |||
47 | /* | ||
48 | Horner evaluation of P=p(A) | ||
49 | d = degree of polynomial; coeff = coefficients of polynomial; W = workspace | ||
50 | */ | ||
51 | 15505 | static PetscErrorCode EvaluatePoly(Mat A,Mat P,Mat W,PetscInt d,PetscScalar *coeff) | |
52 | { | ||
53 | 15505 | PetscInt j; | |
54 | |||
55 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
15505 | PetscFunctionBegin; |
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.
|
15505 | PetscCall(MatZeroEntries(P)); |
57 |
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.
|
15505 | if (!d) PetscCall(MatShift(P,1.0)); |
58 | else { | ||
59 |
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.
|
15505 | PetscCall(MatShift(P,coeff[0])); |
60 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
26264 | for (j=1;j<d;j++) { |
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.
|
10759 | PetscCall(MatMatMult(P,A,MAT_REUSE_MATRIX,PETSC_DEFAULT,&W)); |
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.
|
10759 | PetscCall(MatCopy(W,P,SAME_NONZERO_PATTERN)); |
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.
|
10759 | PetscCall(MatShift(P,coeff[j])); |
64 | } | ||
65 | } | ||
66 |
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.
|
3058 | PetscFunctionReturn(PETSC_SUCCESS); |
67 | } | ||
68 | |||
69 | 12667 | static PetscErrorCode FNEvaluateFunctionMat_Rational(FN fn,Mat A,Mat B) | |
70 | { | ||
71 | 12667 | FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data; | |
72 | 12667 | Mat P,Q,W,F; | |
73 | 12667 | PetscBool iscuda; | |
74 | |||
75 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
12667 | PetscFunctionBegin; |
76 |
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.
|
12667 | if (A==B) PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&P)); |
77 | 12639 | else P = B; | |
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.
|
12667 | PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&W)); |
79 | |||
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.
|
12667 | PetscCall(EvaluatePoly(A,P,W,ctx->np,ctx->pcoeff)); |
81 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
12667 | if (ctx->nq) { |
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.
|
2554 | PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&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.
|
2554 | PetscCall(EvaluatePoly(A,Q,W,ctx->nq,ctx->qcoeff)); |
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.
|
2554 | PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQDENSECUDA,&iscuda)); |
85 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 2 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.
|
5084 | PetscCall(MatGetFactor(Q,iscuda?MATSOLVERCUDA:MATSOLVERPETSC,MAT_FACTOR_LU,&F)); |
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.
|
2554 | PetscCall(MatLUFactorSymbolic(F,Q,NULL,NULL,NULL)); |
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.
|
2554 | PetscCall(MatLUFactorNumeric(F,Q,NULL)); |
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.
|
2554 | PetscCall(MatMatSolve(F,P,P)); |
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.
|
2554 | PetscCall(MatDestroy(&F)); |
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.
|
2554 | PetscCall(MatDestroy(&Q)); |
91 | } | ||
92 | |||
93 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
12667 | if (A==B) { |
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.
|
28 | PetscCall(MatCopy(P,B,SAME_NONZERO_PATTERN)); |
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.
|
28 | PetscCall(MatDestroy(&P)); |
96 | } | ||
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.
|
12667 | PetscCall(MatDestroy(&W)); |
98 |
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.
|
2500 | PetscFunctionReturn(PETSC_SUCCESS); |
99 | } | ||
100 | |||
101 | 256 | static PetscErrorCode FNEvaluateFunctionMatVec_Rational(FN fn,Mat A,Vec v) | |
102 | { | ||
103 | 256 | FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data; | |
104 | 256 | Mat P,Q,W,F; | |
105 | 256 | Vec b; | |
106 | 256 | PetscBool iscuda; | |
107 | |||
108 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
256 | PetscFunctionBegin; |
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.
|
256 | PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&P)); |
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.
|
256 | PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&W)); |
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.
|
256 | PetscCall(EvaluatePoly(A,P,W,ctx->np,ctx->pcoeff)); |
113 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
256 | if (ctx->nq) { |
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.
|
28 | PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q)); |
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.
|
28 | PetscCall(EvaluatePoly(A,Q,W,ctx->nq,ctx->qcoeff)); |
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.
|
28 | PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQDENSECUDA,&iscuda)); |
117 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 2 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.
|
48 | PetscCall(MatGetFactor(Q,iscuda?MATSOLVERCUDA:MATSOLVERPETSC,MAT_FACTOR_LU,&F)); |
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.
|
28 | PetscCall(MatLUFactorSymbolic(F,Q,NULL,NULL,NULL)); |
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.
|
28 | PetscCall(MatLUFactorNumeric(F,Q,NULL)); |
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.
|
28 | PetscCall(MatCreateVecs(P,&b,NULL)); |
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.
|
28 | PetscCall(MatGetColumnVector(P,b,0)); |
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.
|
28 | PetscCall(MatSolve(F,b,v)); |
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.
|
28 | PetscCall(VecDestroy(&b)); |
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.
|
28 | PetscCall(MatDestroy(&F)); |
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.
|
28 | PetscCall(MatDestroy(&Q)); |
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.
|
228 | } else PetscCall(MatGetColumnVector(P,v,0)); |
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.
|
256 | PetscCall(MatDestroy(&P)); |
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.
|
256 | PetscCall(MatDestroy(&W)); |
130 |
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.
|
48 | PetscFunctionReturn(PETSC_SUCCESS); |
131 | } | ||
132 | |||
133 | 33870 | static PetscErrorCode FNEvaluateDerivative_Rational(FN fn,PetscScalar x,PetscScalar *yp) | |
134 | { | ||
135 | 33870 | FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data; | |
136 | 33870 | PetscInt i; | |
137 | 33870 | PetscScalar p,q,pp,qp; | |
138 | |||
139 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
33870 | PetscFunctionBegin; |
140 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
33870 | if (!ctx->np) { |
141 | p = 1.0; | ||
142 | pp = 0.0; | ||
143 | } else { | ||
144 | 33854 | p = ctx->pcoeff[0]; | |
145 | 33854 | pp = 0.0; | |
146 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
54432 | for (i=1;i<ctx->np;i++) { |
147 | 20578 | pp = p+x*pp; | |
148 | 20578 | p = ctx->pcoeff[i]+x*p; | |
149 | } | ||
150 | } | ||
151 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
33870 | if (!ctx->nq) *yp = pp; |
152 | else { | ||
153 | 4423 | q = ctx->qcoeff[0]; | |
154 | 4423 | qp = 0.0; | |
155 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
8902 | for (i=1;i<ctx->nq;i++) { |
156 | 4479 | qp = q+x*qp; | |
157 | 4479 | q = ctx->qcoeff[i]+x*q; | |
158 | } | ||
159 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
4423 | PetscCheck(q!=0.0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Derivative not defined in the requested value"); |
160 | 4423 | *yp = (pp*q-p*qp)/(q*q); | |
161 | } | ||
162 |
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.
|
6607 | PetscFunctionReturn(PETSC_SUCCESS); |
163 | } | ||
164 | |||
165 | 240 | static PetscErrorCode FNView_Rational(FN fn,PetscViewer viewer) | |
166 | { | ||
167 | 240 | FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data; | |
168 | 240 | PetscBool isascii; | |
169 | 240 | PetscInt i; | |
170 | 240 | char str[50]; | |
171 | |||
172 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
240 | PetscFunctionBegin; |
173 |
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.
|
240 | PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); |
174 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
240 | if (isascii) { |
175 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
|
240 | if (fn->alpha!=(PetscScalar)1.0 || fn->beta!=(PetscScalar)1.0) { |
176 |
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.
|
32 | PetscCall(SlepcSNPrintfScalar(str,sizeof(str),fn->alpha,PETSC_FALSE)); |
177 |
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.
|
32 | PetscCall(PetscViewerASCIIPrintf(viewer," scale factors: alpha=%s,",str)); |
178 |
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.
|
32 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); |
179 |
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.
|
32 | PetscCall(SlepcSNPrintfScalar(str,sizeof(str),fn->beta,PETSC_FALSE)); |
180 |
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.
|
32 | PetscCall(PetscViewerASCIIPrintf(viewer," beta=%s\n",str)); |
181 |
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.
|
32 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); |
182 | } | ||
183 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
240 | if (!ctx->nq) { |
184 |
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.
|
150 | if (!ctx->np) PetscCall(PetscViewerASCIIPrintf(viewer," constant: 1.0\n")); |
185 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
150 | else if (ctx->np==1) { |
186 |
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.
|
55 | PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[0],PETSC_FALSE)); |
187 |
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.
|
55 | PetscCall(PetscViewerASCIIPrintf(viewer," constant: %s\n",str)); |
188 | } else { | ||
189 |
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.
|
95 | PetscCall(PetscViewerASCIIPrintf(viewer," polynomial: ")); |
190 |
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.
|
95 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); |
191 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
262 | for (i=0;i<ctx->np-1;i++) { |
192 |
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.
|
167 | PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[i],PETSC_TRUE)); |
193 |
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.
|
167 | PetscCall(PetscViewerASCIIPrintf(viewer,"%s*x^%1" PetscInt_FMT,str,ctx->np-i-1)); |
194 | } | ||
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.
|
95 | PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[ctx->np-1],PETSC_TRUE)); |
196 |
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.
|
95 | PetscCall(PetscViewerASCIIPrintf(viewer,"%s\n",str)); |
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.
|
95 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); |
198 | } | ||
199 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
90 | } else if (!ctx->np) { |
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.
|
16 | PetscCall(PetscViewerASCIIPrintf(viewer," inverse polynomial: 1 / (")); |
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.
|
16 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); |
202 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
48 | for (i=0;i<ctx->nq-1;i++) { |
203 |
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.
|
32 | PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->qcoeff[i],PETSC_TRUE)); |
204 |
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.
|
32 | PetscCall(PetscViewerASCIIPrintf(viewer,"%s*x^%1" PetscInt_FMT,str,ctx->nq-i-1)); |
205 | } | ||
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.
|
16 | PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->qcoeff[ctx->nq-1],PETSC_TRUE)); |
207 |
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.
|
16 | PetscCall(PetscViewerASCIIPrintf(viewer,"%s)\n",str)); |
208 |
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.
|
16 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); |
209 | } else { | ||
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.
|
74 | PetscCall(PetscViewerASCIIPrintf(viewer," rational function: (")); |
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.
|
74 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); |
212 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
148 | for (i=0;i<ctx->np-1;i++) { |
213 |
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.
|
74 | PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[i],PETSC_TRUE)); |
214 |
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.
|
74 | PetscCall(PetscViewerASCIIPrintf(viewer,"%s*x^%1" PetscInt_FMT,str,ctx->np-i-1)); |
215 | } | ||
216 |
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.
|
74 | PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[ctx->np-1],PETSC_TRUE)); |
217 |
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.
|
74 | PetscCall(PetscViewerASCIIPrintf(viewer,"%s) / (",str)); |
218 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
222 | for (i=0;i<ctx->nq-1;i++) { |
219 |
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.
|
148 | PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->qcoeff[i],PETSC_TRUE)); |
220 |
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.
|
148 | PetscCall(PetscViewerASCIIPrintf(viewer,"%s*x^%1" PetscInt_FMT,str,ctx->nq-i-1)); |
221 | } | ||
222 |
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.
|
74 | PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->qcoeff[ctx->nq-1],PETSC_TRUE)); |
223 |
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.
|
74 | PetscCall(PetscViewerASCIIPrintf(viewer,"%s)\n",str)); |
224 |
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.
|
74 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); |
225 | } | ||
226 | } | ||
227 |
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.
|
38 | PetscFunctionReturn(PETSC_SUCCESS); |
228 | } | ||
229 | |||
230 | 2255 | static PetscErrorCode FNRationalSetNumerator_Rational(FN fn,PetscInt np,PetscScalar *pcoeff) | |
231 | { | ||
232 | 2255 | FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data; | |
233 | 2255 | PetscInt i; | |
234 | |||
235 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
2255 | PetscFunctionBegin; |
236 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
2255 | PetscCheck(np>=0,PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Argument np cannot be negative"); |
237 | 2255 | ctx->np = np; | |
238 |
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.
|
2255 | PetscCall(PetscFree(ctx->pcoeff)); |
239 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2255 | if (np) { |
240 |
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.
|
2239 | PetscCall(PetscMalloc1(np,&ctx->pcoeff)); |
241 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
5872 | for (i=0;i<np;i++) ctx->pcoeff[i] = pcoeff[i]; |
242 | } | ||
243 |
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.
|
430 | PetscFunctionReturn(PETSC_SUCCESS); |
244 | } | ||
245 | |||
246 | /*@ | ||
247 | FNRationalSetNumerator - Sets the parameters defining the numerator of the | ||
248 | rational function. | ||
249 | |||
250 | Logically Collective | ||
251 | |||
252 | Input Parameters: | ||
253 | + fn - the math function context | ||
254 | . np - number of coefficients | ||
255 | - pcoeff - coefficients (array of scalar values) | ||
256 | |||
257 | Notes: | ||
258 | Let the rational function r(x) = p(x)/q(x), where p(x) and q(x) are polynomials. | ||
259 | This function provides the coefficients of the numerator p(x). | ||
260 | Hence, p(x) is of degree np-1. | ||
261 | If np is zero, then the numerator is assumed to be p(x)=1. | ||
262 | |||
263 | In polynomials, high order coefficients are stored in the first positions | ||
264 | of the array, e.g. to represent x^2-3 use {1,0,-3}. | ||
265 | |||
266 | Level: intermediate | ||
267 | |||
268 | .seealso: FNRationalSetDenominator(), FNRationalGetNumerator() | ||
269 | @*/ | ||
270 | 2255 | PetscErrorCode FNRationalSetNumerator(FN fn,PetscInt np,PetscScalar pcoeff[]) | |
271 | { | ||
272 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
2255 | PetscFunctionBegin; |
273 |
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.
|
2255 | PetscValidHeaderSpecific(fn,FN_CLASSID,1); |
274 |
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.
|
2255 | PetscValidLogicalCollectiveInt(fn,np,2); |
275 |
4/10✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
2255 | if (np) PetscAssertPointer(pcoeff,3); |
276 |
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.
|
2255 | PetscTryMethod(fn,"FNRationalSetNumerator_C",(FN,PetscInt,PetscScalar*),(fn,np,pcoeff)); |
277 |
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.
|
2255 | PetscFunctionReturn(PETSC_SUCCESS); |
278 | } | ||
279 | |||
280 | 16 | static PetscErrorCode FNRationalGetNumerator_Rational(FN fn,PetscInt *np,PetscScalar *pcoeff[]) | |
281 | { | ||
282 | 16 | FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data; | |
283 | 16 | PetscInt i; | |
284 | |||
285 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
16 | PetscFunctionBegin; |
286 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
16 | if (np) *np = ctx->np; |
287 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
16 | if (pcoeff) { |
288 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
16 | if (!ctx->np) *pcoeff = NULL; |
289 | else { | ||
290 |
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.
|
16 | PetscCall(PetscMalloc1(ctx->np,pcoeff)); |
291 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
48 | for (i=0;i<ctx->np;i++) (*pcoeff)[i] = ctx->pcoeff[i]; |
292 | } | ||
293 | } | ||
294 |
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.
|
2 | PetscFunctionReturn(PETSC_SUCCESS); |
295 | } | ||
296 | |||
297 | /*@C | ||
298 | FNRationalGetNumerator - Gets the parameters that define the numerator of the | ||
299 | rational function. | ||
300 | |||
301 | Not Collective | ||
302 | |||
303 | Input Parameter: | ||
304 | . fn - the math function context | ||
305 | |||
306 | Output Parameters: | ||
307 | + np - number of coefficients | ||
308 | - pcoeff - coefficients (array of scalar values, length nq) | ||
309 | |||
310 | Notes: | ||
311 | The values passed by user with FNRationalSetNumerator() are returned (or null | ||
312 | pointers otherwise). | ||
313 | The pcoeff array should be freed by the user when no longer needed. | ||
314 | |||
315 | Level: intermediate | ||
316 | |||
317 | .seealso: FNRationalSetNumerator() | ||
318 | @*/ | ||
319 | 16 | PetscErrorCode FNRationalGetNumerator(FN fn,PetscInt *np,PetscScalar *pcoeff[]) PeNS | |
320 | { | ||
321 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
16 | PetscFunctionBegin; |
322 |
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.
|
16 | PetscValidHeaderSpecific(fn,FN_CLASSID,1); |
323 |
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.
|
16 | PetscUseMethod(fn,"FNRationalGetNumerator_C",(FN,PetscInt*,PetscScalar**),(fn,np,pcoeff)); |
324 |
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.
|
16 | PetscFunctionReturn(PETSC_SUCCESS); |
325 | } | ||
326 | |||
327 | 352 | static PetscErrorCode FNRationalSetDenominator_Rational(FN fn,PetscInt nq,PetscScalar *qcoeff) | |
328 | { | ||
329 | 352 | FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data; | |
330 | 352 | PetscInt i; | |
331 | |||
332 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
352 | PetscFunctionBegin; |
333 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
352 | PetscCheck(nq>=0,PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Argument nq cannot be negative"); |
334 | 352 | ctx->nq = nq; | |
335 |
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.
|
352 | PetscCall(PetscFree(ctx->qcoeff)); |
336 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
352 | if (nq) { |
337 |
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(PetscMalloc1(nq,&ctx->qcoeff)); |
338 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1098 | for (i=0;i<nq;i++) ctx->qcoeff[i] = qcoeff[i]; |
339 | } | ||
340 |
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.
|
65 | PetscFunctionReturn(PETSC_SUCCESS); |
341 | } | ||
342 | |||
343 | /*@ | ||
344 | FNRationalSetDenominator - Sets the parameters defining the denominator of the | ||
345 | rational function. | ||
346 | |||
347 | Logically Collective | ||
348 | |||
349 | Input Parameters: | ||
350 | + fn - the math function context | ||
351 | . nq - number of coefficients | ||
352 | - qcoeff - coefficients (array of scalar values) | ||
353 | |||
354 | Notes: | ||
355 | Let the rational function r(x) = p(x)/q(x), where p(x) and q(x) are polynomials. | ||
356 | This function provides the coefficients of the denominator q(x). | ||
357 | Hence, q(x) is of degree nq-1. | ||
358 | If nq is zero, then the function is assumed to be polynomial, r(x) = p(x). | ||
359 | |||
360 | In polynomials, high order coefficients are stored in the first positions | ||
361 | of the array, e.g. to represent x^2-3 use {1,0,-3}. | ||
362 | |||
363 | Level: intermediate | ||
364 | |||
365 | .seealso: FNRationalSetNumerator(), FNRationalGetDenominator() | ||
366 | @*/ | ||
367 | 352 | PetscErrorCode FNRationalSetDenominator(FN fn,PetscInt nq,PetscScalar qcoeff[]) | |
368 | { | ||
369 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
352 | PetscFunctionBegin; |
370 |
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.
|
352 | PetscValidHeaderSpecific(fn,FN_CLASSID,1); |
371 |
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.
|
352 | PetscValidLogicalCollectiveInt(fn,nq,2); |
372 |
4/10✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
352 | if (nq) PetscAssertPointer(qcoeff,3); |
373 |
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.
|
352 | PetscTryMethod(fn,"FNRationalSetDenominator_C",(FN,PetscInt,PetscScalar*),(fn,nq,qcoeff)); |
374 |
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.
|
352 | PetscFunctionReturn(PETSC_SUCCESS); |
375 | } | ||
376 | |||
377 | 116 | static PetscErrorCode FNRationalGetDenominator_Rational(FN fn,PetscInt *nq,PetscScalar *qcoeff[]) | |
378 | { | ||
379 | 116 | FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data; | |
380 | 116 | PetscInt i; | |
381 | |||
382 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
116 | PetscFunctionBegin; |
383 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
116 | if (nq) *nq = ctx->nq; |
384 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
116 | if (qcoeff) { |
385 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
116 | if (!ctx->nq) *qcoeff = NULL; |
386 | else { | ||
387 |
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.
|
46 | PetscCall(PetscMalloc1(ctx->nq,qcoeff)); |
388 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
154 | for (i=0;i<ctx->nq;i++) (*qcoeff)[i] = ctx->qcoeff[i]; |
389 | } | ||
390 | } | ||
391 |
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.
|
22 | PetscFunctionReturn(PETSC_SUCCESS); |
392 | } | ||
393 | |||
394 | /*@C | ||
395 | FNRationalGetDenominator - Gets the parameters that define the denominator of the | ||
396 | rational function. | ||
397 | |||
398 | Not Collective | ||
399 | |||
400 | Input Parameter: | ||
401 | . fn - the math function context | ||
402 | |||
403 | Output Parameters: | ||
404 | + nq - number of coefficients | ||
405 | - qcoeff - coefficients (array of scalar values, length nq) | ||
406 | |||
407 | Notes: | ||
408 | The values passed by user with FNRationalSetDenominator() are returned (or a null | ||
409 | pointer otherwise). | ||
410 | The qcoeff array should be freed by the user when no longer needed. | ||
411 | |||
412 | Level: intermediate | ||
413 | |||
414 | .seealso: FNRationalSetDenominator() | ||
415 | @*/ | ||
416 | 116 | PetscErrorCode FNRationalGetDenominator(FN fn,PetscInt *nq,PetscScalar *qcoeff[]) PeNS | |
417 | { | ||
418 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
116 | PetscFunctionBegin; |
419 |
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.
|
116 | PetscValidHeaderSpecific(fn,FN_CLASSID,1); |
420 |
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.
|
116 | PetscUseMethod(fn,"FNRationalGetDenominator_C",(FN,PetscInt*,PetscScalar**),(fn,nq,qcoeff)); |
421 |
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.
|
116 | PetscFunctionReturn(PETSC_SUCCESS); |
422 | } | ||
423 | |||
424 | 2107 | static PetscErrorCode FNSetFromOptions_Rational(FN fn,PetscOptionItems PetscOptionsObject) | |
425 | { | ||
426 | #define PARMAX 10 | ||
427 | 2107 | PetscScalar array[PARMAX]; | |
428 | 2107 | PetscInt i,k; | |
429 | 2107 | PetscBool flg; | |
430 | |||
431 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
2107 | PetscFunctionBegin; |
432 |
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.
|
2107 | PetscOptionsHeadBegin(PetscOptionsObject,"FN Rational Options"); |
433 | |||
434 | 2107 | k = PARMAX; | |
435 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
23177 | for (i=0;i<k;i++) array[i] = 0; |
436 |
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.
|
2107 | PetscCall(PetscOptionsScalarArray("-fn_rational_numerator","Numerator coefficients (one or more scalar values separated with a comma without spaces)","FNRationalSetNumerator",array,&k,&flg)); |
437 |
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.
|
2107 | if (flg) PetscCall(FNRationalSetNumerator(fn,k,array)); |
438 | |||
439 | 2107 | k = PARMAX; | |
440 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
23177 | for (i=0;i<k;i++) array[i] = 0; |
441 |
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.
|
2107 | PetscCall(PetscOptionsScalarArray("-fn_rational_denominator","Denominator coefficients (one or more scalar values separated with a comma without spaces)","FNRationalSetDenominator",array,&k,&flg)); |
442 |
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.
|
2107 | if (flg) PetscCall(FNRationalSetDenominator(fn,k,array)); |
443 | |||
444 |
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.
|
2107 | PetscOptionsHeadEnd(); |
445 |
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.
|
402 | PetscFunctionReturn(PETSC_SUCCESS); |
446 | } | ||
447 | |||
448 | 288 | static PetscErrorCode FNDuplicate_Rational(FN fn,MPI_Comm comm,FN *newfn) | |
449 | { | ||
450 | 288 | FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data,*ctx2 = (FN_RATIONAL*)(*newfn)->data; | |
451 | 288 | PetscInt i; | |
452 | |||
453 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
288 | PetscFunctionBegin; |
454 | 288 | ctx2->np = ctx->np; | |
455 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
288 | if (ctx->np) { |
456 |
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(PetscMalloc1(ctx->np,&ctx2->pcoeff)); |
457 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
808 | for (i=0;i<ctx->np;i++) ctx2->pcoeff[i] = ctx->pcoeff[i]; |
458 | } | ||
459 | 288 | ctx2->nq = ctx->nq; | |
460 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
288 | if (ctx->nq) { |
461 |
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.
|
104 | PetscCall(PetscMalloc1(ctx->nq,&ctx2->qcoeff)); |
462 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
336 | for (i=0;i<ctx->nq;i++) ctx2->qcoeff[i] = ctx->qcoeff[i]; |
463 | } | ||
464 |
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.
|
56 | PetscFunctionReturn(PETSC_SUCCESS); |
465 | } | ||
466 | |||
467 | 2519 | static PetscErrorCode FNDestroy_Rational(FN fn) | |
468 | { | ||
469 | 2519 | FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data; | |
470 | |||
471 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
2519 | PetscFunctionBegin; |
472 |
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.
|
2519 | PetscCall(PetscFree(ctx->pcoeff)); |
473 |
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.
|
2519 | PetscCall(PetscFree(ctx->qcoeff)); |
474 |
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.
|
2519 | PetscCall(PetscFree(fn->data)); |
475 |
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.
|
2519 | PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetNumerator_C",NULL)); |
476 |
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.
|
2519 | PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetNumerator_C",NULL)); |
477 |
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.
|
2519 | PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetDenominator_C",NULL)); |
478 |
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.
|
2519 | PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetDenominator_C",NULL)); |
479 |
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.
|
482 | PetscFunctionReturn(PETSC_SUCCESS); |
480 | } | ||
481 | |||
482 | 2519 | SLEPC_EXTERN PetscErrorCode FNCreate_Rational(FN fn) | |
483 | { | ||
484 | 2519 | FN_RATIONAL *ctx; | |
485 | |||
486 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
2519 | PetscFunctionBegin; |
487 |
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.
|
2519 | PetscCall(PetscNew(&ctx)); |
488 | 2519 | fn->data = (void*)ctx; | |
489 | |||
490 | 2519 | fn->ops->evaluatefunction = FNEvaluateFunction_Rational; | |
491 | 2519 | fn->ops->evaluatederivative = FNEvaluateDerivative_Rational; | |
492 | 2519 | fn->ops->evaluatefunctionmat[0] = FNEvaluateFunctionMat_Rational; | |
493 | 2519 | fn->ops->evaluatefunctionmatvec[0] = FNEvaluateFunctionMatVec_Rational; | |
494 | #if defined(PETSC_HAVE_CUDA) | ||
495 | 525 | fn->ops->evaluatefunctionmatcuda[0] = FNEvaluateFunctionMat_Rational; | |
496 | 525 | fn->ops->evaluatefunctionmatveccuda[0] = FNEvaluateFunctionMatVec_Rational; | |
497 | #endif | ||
498 | 2519 | fn->ops->setfromoptions = FNSetFromOptions_Rational; | |
499 | 2519 | fn->ops->view = FNView_Rational; | |
500 | 2519 | fn->ops->duplicate = FNDuplicate_Rational; | |
501 | 2519 | fn->ops->destroy = FNDestroy_Rational; | |
502 |
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.
|
2519 | PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetNumerator_C",FNRationalSetNumerator_Rational)); |
503 |
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.
|
2519 | PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetNumerator_C",FNRationalGetNumerator_Rational)); |
504 |
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.
|
2519 | PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetDenominator_C",FNRationalSetDenominator_Rational)); |
505 |
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.
|
2519 | PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetDenominator_C",FNRationalGetDenominator_Rational)); |
506 |
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.
|
482 | PetscFunctionReturn(PETSC_SUCCESS); |
507 | } | ||
508 |