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