1 /* 2 Common tools for constructing discretizations 3 */ 4 #if !defined(PETSCDT_H) 5 #define PETSCDT_H 6 7 #include <petscsys.h> 8 9 PETSC_EXTERN PetscClassId PETSCQUADRATURE_CLASSID; 10 11 /*S 12 PetscQuadrature - Quadrature rule for integration. 13 14 Level: beginner 15 16 .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy() 17 S*/ 18 typedef struct _p_PetscQuadrature *PetscQuadrature; 19 20 /*E 21 PetscGaussLobattoLegendreCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights 22 23 Level: intermediate 24 25 $ PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA - compute the nodes via linear algebra 26 $ PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON - compute the nodes by solving a nonlinear equation with Newton's method 27 28 E*/ 29 typedef enum {PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON} PetscGaussLobattoLegendreCreateType; 30 31 PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *); 32 PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *); 33 PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt*); 34 PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt); 35 PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt*); 36 PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt); 37 PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt*, PetscInt*, PetscInt*, const PetscReal *[], const PetscReal *[]); 38 PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal [], const PetscReal []); 39 PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer); 40 PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *); 41 42 PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *); 43 44 PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *); 45 46 PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*); 47 PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt,PetscReal,PetscReal,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*); 48 PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt,PetscReal,PetscReal,PetscReal*,PetscReal*); 49 PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*,PetscReal*); 50 PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*,PetscReal*); 51 PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt,PetscGaussLobattoLegendreCreateType,PetscReal*,PetscReal*); 52 PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt,PetscInt,const PetscReal*,PetscInt,const PetscReal*,PetscReal*); 53 PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*); 54 PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*); 55 56 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *); 57 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *); 58 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *); 59 60 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *); 61 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 62 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 63 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 64 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 65 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 66 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 67 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 68 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 69 70 PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 71 PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 72 PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 73 PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *); 74 PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *); 75 PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 76 PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *); 77 PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]); 78 PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 79 80 #if defined(PETSC_USE_64BIT_INDICES) 81 #define PETSC_FACTORIAL_MAX 20 82 #define PETSC_BINOMIAL_MAX 61 83 #else 84 #define PETSC_FACTORIAL_MAX 12 85 #define PETSC_BINOMIAL_MAX 29 86 #endif 87 88 /*MC 89 PetscDTFactorial - Approximate n! as a real number 90 91 Input Arguments: 92 . n - a non-negative integer 93 94 Output Arguments: 95 . factorial - n! 96 97 Level: beginner 98 M*/ 99 PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial) 100 { 101 PetscReal f = 1.0; 102 PetscInt i; 103 104 PetscFunctionBegin; 105 *factorial = -1.0; 106 if (n < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %D\n", n); 107 for (i = 1; i < n+1; ++i) f *= (PetscReal)i; 108 *factorial = f; 109 PetscFunctionReturn(0); 110 } 111 112 /*MC 113 PetscDTFactorialInt - Compute n! as an integer 114 115 Input Arguments: 116 . n - a non-negative integer 117 118 Output Arguments: 119 . factorial - n! 120 121 Level: beginner 122 123 Note: this is limited to n such that n! can be represented by PetscInt, which is 12 if PetscInt is a signed 32-bit integer and 20 if PetscInt is a signed 64-bit integer. 124 M*/ 125 PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial) 126 { 127 PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600}; 128 129 PetscFunctionBegin; 130 *factorial = -1; 131 if (n < 0 || n > PETSC_FACTORIAL_MAX) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of elements %D is not in supported range [0,%D]\n",n,PETSC_FACTORIAL_MAX); 132 if (n <= 12) { 133 *factorial = facLookup[n]; 134 } else { 135 PetscInt f = facLookup[12]; 136 PetscInt i; 137 138 for (i = 13; i < n+1; ++i) f *= i; 139 *factorial = f; 140 } 141 PetscFunctionReturn(0); 142 } 143 144 /*MC 145 PetscDTBinomial - Approximate the binomial coefficient "n choose k" 146 147 Input Arguments: 148 + n - a non-negative integer 149 - k - an integer between 0 and n, inclusive 150 151 Output Arguments: 152 . binomial - approximation of the binomial coefficient n choose k 153 154 Level: beginner 155 M*/ 156 PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial) 157 { 158 PetscFunctionBeginHot; 159 *binomial = -1.0; 160 if (n < 0 || k < 0 || k > n) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial arguments (%D %D) must be non-negative, k <= n\n", n, k); 161 if (n <= 3) { 162 PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}}; 163 164 *binomial = (PetscReal)binomLookup[n][k]; 165 } else { 166 PetscReal binom = 1.0; 167 PetscInt i; 168 169 k = PetscMin(k, n - k); 170 for (i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1); 171 *binomial = binom; 172 } 173 PetscFunctionReturn(0); 174 } 175 176 /*MC 177 PetscDTBinomialInt - Compute the binomial coefficient "n choose k" 178 179 Input Arguments: 180 + n - a non-negative integer 181 - k - an integer between 0 and n, inclusive 182 183 Output Arguments: 184 . binomial - the binomial coefficient n choose k 185 186 Note: this is limited by integers that can be represented by PetscInt 187 188 Level: beginner 189 M*/ 190 PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial) 191 { 192 PetscInt bin; 193 194 PetscFunctionBegin; 195 *binomial = -1; 196 if (n < 0 || k < 0 || k > n) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial arguments (%D %D) must be non-negative, k <= n\n", n, k); 197 if (n > PETSC_BINOMIAL_MAX) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial elements %D is larger than max for PetscInt, %D\n", n, PETSC_BINOMIAL_MAX); 198 if (n <= 3) { 199 PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}}; 200 201 bin = binomLookup[n][k]; 202 } else { 203 PetscInt binom = 1; 204 PetscInt i; 205 206 k = PetscMin(k, n - k); 207 for (i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1); 208 bin = binom; 209 } 210 *binomial = bin; 211 PetscFunctionReturn(0); 212 } 213 214 /*MC 215 PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps. 216 217 A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation, 218 by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in 219 some position j >= i. This swap is encoded as the difference (j - i). The difference d_i at step i is less than 220 (n - i). This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number 221 (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}. 222 223 Input Arguments: 224 + n - a non-negative integer (see note about limits below) 225 - k - an integer in [0, n!) 226 227 Output Arguments: 228 + perm - the permuted list of the integers [0, ..., n-1] 229 - isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps. 230 231 Note: this is limited to n such that n! can be represented by PetscInt, which is 12 if PetscInt is a signed 32-bit integer and 20 if PetscInt is a signed 64-bit integer. 232 233 Level: beginner 234 M*/ 235 PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd) 236 { 237 PetscInt odd = 0; 238 PetscInt i; 239 PetscInt work[PETSC_FACTORIAL_MAX]; 240 PetscInt *w; 241 242 PetscFunctionBegin; 243 if (isOdd) *isOdd = PETSC_FALSE; 244 if (n < 0 || n > PETSC_FACTORIAL_MAX) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of elements %D is not in supported range [0,%D]\n",n,PETSC_FACTORIAL_MAX); 245 w = &work[n - 2]; 246 for (i = 2; i <= n; i++) { 247 *(w--) = k % i; 248 k /= i; 249 } 250 for (i = 0; i < n; i++) perm[i] = i; 251 for (i = 0; i < n - 1; i++) { 252 PetscInt s = work[i]; 253 PetscInt swap = perm[i]; 254 255 perm[i] = perm[i + s]; 256 perm[i + s] = swap; 257 odd ^= (!!s); 258 } 259 if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 260 PetscFunctionReturn(0); 261 } 262 263 /*MC 264 PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!). This inverts PetscDTEnumPerm. 265 266 Input Arguments: 267 + n - a non-negative integer (see note about limits below) 268 - perm - the permuted list of the integers [0, ..., n-1] 269 270 Output Arguments: 271 + k - an integer in [0, n!) 272 - isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps. 273 274 Note: this is limited to n such that n! can be represented by PetscInt, which is 12 if PetscInt is a signed 32-bit integer and 20 if PetscInt is a signed 64-bit integer. 275 276 Level: beginner 277 M*/ 278 PETSC_STATIC_INLINE PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd) 279 { 280 PetscInt odd = 0; 281 PetscInt i, idx; 282 PetscInt work[PETSC_FACTORIAL_MAX]; 283 PetscInt iwork[PETSC_FACTORIAL_MAX]; 284 285 PetscFunctionBeginHot; 286 *k = -1; 287 if (isOdd) *isOdd = PETSC_FALSE; 288 if (n < 0 || n > PETSC_FACTORIAL_MAX) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of elements %D is not in supported range [0,%D]\n",n,PETSC_FACTORIAL_MAX); 289 for (i = 0; i < n; i++) work[i] = i; /* partial permutation */ 290 for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */ 291 for (idx = 0, i = 0; i < n - 1; i++) { 292 PetscInt j = perm[i]; 293 PetscInt icur = work[i]; 294 PetscInt jloc = iwork[j]; 295 PetscInt diff = jloc - i; 296 297 idx = idx * (n - i) + diff; 298 /* swap (i, jloc) */ 299 work[i] = j; 300 work[jloc] = icur; 301 iwork[j] = i; 302 iwork[icur] = jloc; 303 odd ^= (!!diff); 304 } 305 *k = idx; 306 if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 307 PetscFunctionReturn(0); 308 } 309 310 /*MC 311 PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k). 312 The encoding is in lexicographic order. 313 314 Input Arguments: 315 + n - a non-negative integer (see note about limits below) 316 . k - an integer in [0, n] 317 - j - an index in [0, n choose k) 318 319 Output Arguments: 320 . subset - the jth subset of size k of the integers [0, ..., n - 1] 321 322 Note: this is limited by arguments such that n choose k can be represented by PetscInt 323 324 Level: beginner 325 326 .seealso: PetscDTSubsetIndex() 327 M*/ 328 PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset) 329 { 330 PetscInt Nk, i, l; 331 PetscErrorCode ierr; 332 333 PetscFunctionBeginHot; 334 ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr); 335 for (i = 0, l = 0; i < n && l < k; i++) { 336 PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 337 PetscInt Nminusk = Nk - Nminuskminus; 338 339 if (j < Nminuskminus) { 340 subset[l++] = i; 341 Nk = Nminuskminus; 342 } else { 343 j -= Nminuskminus; 344 Nk = Nminusk; 345 } 346 } 347 PetscFunctionReturn(0); 348 } 349 350 /*MC 351 PetscDTSubsetIndex - Convert an ordered subset of k integers from the set [0, ..., n - 1] to its encoding as an integers in [0, n choose k) in lexicographic order. This is the inverse of PetscDTEnumSubset. 352 353 Input Arguments: 354 + n - a non-negative integer (see note about limits below) 355 . k - an integer in [0, n] 356 - subset - an ordered subset of the integers [0, ..., n - 1] 357 358 Output Arguments: 359 . index - the rank of the subset in lexicographic order 360 361 Note: this is limited by arguments such that n choose k can be represented by PetscInt 362 363 Level: beginner 364 365 .seealso: PetscDTEnumSubset() 366 M*/ 367 PETSC_STATIC_INLINE PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index) 368 { 369 PetscInt i, j = 0, l, Nk; 370 PetscErrorCode ierr; 371 372 PetscFunctionBegin; 373 *index = -1; 374 ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr); 375 for (i = 0, l = 0; i < n && l < k; i++) { 376 PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 377 PetscInt Nminusk = Nk - Nminuskminus; 378 379 if (subset[l] == i) { 380 l++; 381 Nk = Nminuskminus; 382 } else { 383 j += Nminuskminus; 384 Nk = Nminusk; 385 } 386 } 387 *index = j; 388 PetscFunctionReturn(0); 389 } 390 391 /*MC 392 PetscDTEnumSubset - Split the integers [0, ..., n - 1] into two complementary ordered subsets, the first subset of size k and being the jth subset of that size in lexicographic order. 393 394 Input Arguments: 395 + n - a non-negative integer (see note about limits below) 396 . k - an integer in [0, n] 397 - j - an index in [0, n choose k) 398 399 Output Arguments: 400 + perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set. 401 - isOdd - if not NULL, return whether perm is an even or odd permutation. 402 403 Note: this is limited by arguments such that n choose k can be represented by PetscInt 404 405 Level: beginner 406 407 .seealso: PetscDTEnumSubset(), PetscDTSubsetIndex() 408 M*/ 409 PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd) 410 { 411 PetscInt i, l, m, *subcomp, Nk; 412 PetscInt odd; 413 PetscErrorCode ierr; 414 415 PetscFunctionBegin; 416 if (isOdd) *isOdd = PETSC_FALSE; 417 ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr); 418 odd = 0; 419 subcomp = &perm[k]; 420 for (i = 0, l = 0, m = 0; i < n && l < k; i++) { 421 PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 422 PetscInt Nminusk = Nk - Nminuskminus; 423 424 if (j < Nminuskminus) { 425 perm[l++] = i; 426 Nk = Nminuskminus; 427 } else { 428 subcomp[m++] = i; 429 j -= Nminuskminus; 430 odd ^= ((k - l) & 1); 431 Nk = Nminusk; 432 } 433 } 434 for (; i < n; i++) { 435 subcomp[m++] = i; 436 } 437 if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 438 PetscFunctionReturn(0); 439 } 440 441 struct _p_PetscTabulation { 442 PetscInt K; /* Indicates a k-jet, namely tabulated derviatives up to order k */ 443 PetscInt Nr; /* THe number of tabulation replicas (often 1) */ 444 PetscInt Np; /* The number of tabulation points in a replica */ 445 PetscInt Nb; /* The number of functions tabulated */ 446 PetscInt Nc; /* The number of function components */ 447 PetscInt cdim; /* The coordinate dimension */ 448 PetscReal **T; /* The tabulation T[K] of functions and their derivatives 449 T[0] = B[Nr*Np][Nb][Nc]: The basis function values at quadrature points 450 T[1] = D[Nr*Np][Nb][Nc][cdim]: The basis function derivatives at quadrature points 451 T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */ 452 }; 453 typedef struct _p_PetscTabulation *PetscTabulation; 454 455 #endif 456