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