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