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