1 /* 2 Common tools for constructing discretizations 3 */ 4 #pragma once 5 6 #include <petscsys.h> 7 #include <petscdmtypes.h> 8 #include <petscistypes.h> 9 10 /* MANSEC = DM */ 11 /* SUBMANSEC = DT */ 12 13 PETSC_EXTERN PetscClassId PETSCQUADRATURE_CLASSID; 14 15 /*S 16 PetscQuadrature - Quadrature rule for numerical integration. 17 18 Level: beginner 19 20 .seealso: `PetscQuadratureCreate()`, `PetscQuadratureDestroy()` 21 S*/ 22 typedef struct _p_PetscQuadrature *PetscQuadrature; 23 24 /*E 25 PetscGaussLobattoLegendreCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights 26 27 Values: 28 + `PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA` - compute the nodes via linear algebra 29 - `PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON` - compute the nodes by solving a nonlinear equation with Newton's method 30 31 Level: intermediate 32 33 .seealso: `PetscQuadrature` 34 E*/ 35 typedef enum { 36 PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA, 37 PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON 38 } PetscGaussLobattoLegendreCreateType; 39 40 /*E 41 PetscDTNodeType - A description of strategies for generating nodes (both 42 quadrature nodes and nodes for Lagrange polynomials) 43 44 Values: 45 + `PETSCDTNODES_DEFAULT` - Nodes chosen by PETSc 46 . `PETSCDTNODES_GAUSSJACOBI` - Nodes at either Gauss-Jacobi or Gauss-Lobatto-Jacobi quadrature points 47 . `PETSCDTNODES_EQUISPACED` - Nodes equispaced either including the endpoints or excluding them 48 - `PETSCDTNODES_TANHSINH` - Nodes at Tanh-Sinh quadrature points 49 50 Level: intermediate 51 52 Note: 53 A `PetscDTNodeType` can be paired with a `PetscBool` to indicate whether 54 the nodes include endpoints or not, and in the case of `PETSCDT_GAUSSJACOBI` 55 with exponents for the weight function. 56 57 .seealso: `PetscQuadrature` 58 E*/ 59 typedef enum { 60 PETSCDTNODES_DEFAULT = -1, 61 PETSCDTNODES_GAUSSJACOBI = 0, 62 PETSCDTNODES_EQUISPACED = 1, 63 PETSCDTNODES_TANHSINH = 2 64 } PetscDTNodeType; 65 66 PETSC_EXTERN const char *const *const PetscDTNodeTypes; 67 68 /*E 69 PetscDTSimplexQuadratureType - A description of classes of quadrature rules for simplices 70 71 Values: 72 + `PETSCDTSIMPLEXQUAD_DEFAULT` - Quadrature rule chosen by PETSc 73 . `PETSCDTSIMPLEXQUAD_CONIC` - Quadrature rules constructed as 74 conically-warped tensor products of 1D 75 Gauss-Jacobi quadrature rules. These are 76 explicitly computable in any dimension for any 77 degree, and the tensor-product structure can be 78 exploited by sum-factorization methods, but 79 they are not efficient in terms of nodes per 80 polynomial degree. 81 - `PETSCDTSIMPLEXQUAD_MINSYM` - Quadrature rules that are fully symmetric 82 (symmetries of the simplex preserve the nodes 83 and weights) with minimal (or near minimal) 84 number of nodes. In dimensions higher than 1 85 these are not simple to compute, so lookup 86 tables are used. 87 88 Level: intermediate 89 90 .seealso: `PetscQuadrature`, `PetscDTSimplexQuadrature()` 91 E*/ 92 typedef enum { 93 PETSCDTSIMPLEXQUAD_DEFAULT = -1, 94 PETSCDTSIMPLEXQUAD_CONIC = 0, 95 PETSCDTSIMPLEXQUAD_MINSYM = 1 96 } PetscDTSimplexQuadratureType; 97 98 PETSC_EXTERN const char *const *const PetscDTSimplexQuadratureTypes; 99 100 PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *); 101 PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *); 102 PETSC_EXTERN PetscErrorCode PetscQuadratureGetCellType(PetscQuadrature, DMPolytopeType *); 103 PETSC_EXTERN PetscErrorCode PetscQuadratureSetCellType(PetscQuadrature, DMPolytopeType); 104 PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt *); 105 PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt); 106 PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt *); 107 PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt); 108 PETSC_EXTERN PetscErrorCode PetscQuadratureEqual(PetscQuadrature, PetscQuadrature, PetscBool *); 109 PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt *, PetscInt *, PetscInt *, const PetscReal *[], const PetscReal *[]); 110 PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal[], const PetscReal[]); 111 PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer); 112 PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *); 113 114 PETSC_EXTERN PetscErrorCode PetscDTTensorQuadratureCreate(PetscQuadrature, PetscQuadrature, PetscQuadrature *); 115 PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *); 116 PETSC_EXTERN PetscErrorCode PetscQuadratureComputePermutations(PetscQuadrature, PetscInt *, IS *[]); 117 118 PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *); 119 120 PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *); 121 PETSC_EXTERN PetscErrorCode PetscDTJacobiNorm(PetscReal, PetscReal, PetscInt, PetscReal *); 122 PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt, PetscReal, PetscReal, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *); 123 PETSC_EXTERN PetscErrorCode PetscDTJacobiEvalJet(PetscReal, PetscReal, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]); 124 PETSC_EXTERN PetscErrorCode PetscDTPKDEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]); 125 PETSC_EXTERN PetscErrorCode PetscDTPTrimmedSize(PetscInt, PetscInt, PetscInt, PetscInt *); 126 PETSC_EXTERN PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscInt, PetscReal[]); 127 PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt, PetscReal, PetscReal, PetscReal *, PetscReal *); 128 PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *); 129 PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *); 130 PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt, PetscGaussLobattoLegendreCreateType, PetscReal *, PetscReal *); 131 PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *); 132 PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *); 133 PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *); 134 PETSC_EXTERN PetscErrorCode PetscDTSimplexQuadrature(PetscInt, PetscInt, PetscDTSimplexQuadratureType, PetscQuadrature *); 135 PETSC_EXTERN PetscErrorCode PetscDTCreateDefaultQuadrature(DMPolytopeType, PetscInt, PetscQuadrature *, PetscQuadrature *); 136 137 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *); 138 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *); 139 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *); 140 141 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *); 142 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 143 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 144 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 145 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 146 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 147 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 148 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 149 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 150 151 /*MC 152 PETSC_FORM_DEGREE_UNDEFINED - Indicates that a field does not have 153 a well-defined form degree in exterior calculus. 154 155 Level: advanced 156 157 .seealso: `PetscDTAltV`, `PetscDualSpaceGetFormDegree()` 158 M*/ 159 #define PETSC_FORM_DEGREE_UNDEFINED PETSC_INT_MIN 160 161 PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 162 PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 163 PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 164 PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *); 165 PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *); 166 PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 167 PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *); 168 PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]); 169 PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 170 171 PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt, PetscInt, const PetscInt[], PetscInt *); 172 PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt, PetscInt, PetscInt, PetscInt[]); 173 PETSC_EXTERN PetscErrorCode PetscDTGradedOrderToIndex(PetscInt, const PetscInt[], PetscInt *); 174 PETSC_EXTERN PetscErrorCode PetscDTIndexToGradedOrder(PetscInt, PetscInt, PetscInt[]); 175 176 #if defined(PETSC_USE_64BIT_INDICES) 177 #define PETSC_FACTORIAL_MAX 20 178 #define PETSC_BINOMIAL_MAX 61 179 #else 180 #define PETSC_FACTORIAL_MAX 12 181 #define PETSC_BINOMIAL_MAX 29 182 #endif 183 184 /*MC 185 PetscDTFactorial - Approximate n! as a real number 186 187 Input Parameter: 188 . n - a non-negative integer 189 190 Output Parameter: 191 . factorial - n! 192 193 Level: beginner 194 195 .seealso: `PetscDTFactorialInt()`, `PetscDTBinomialInt()`, `PetscDTBinomial()` 196 M*/ 197 static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial) 198 { 199 PetscReal f = 1.0; 200 201 PetscFunctionBegin; 202 *factorial = -1.0; 203 PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %" PetscInt_FMT, n); 204 for (PetscInt i = 1; i < n + 1; ++i) f *= (PetscReal)i; 205 *factorial = f; 206 PetscFunctionReturn(PETSC_SUCCESS); 207 } 208 209 /*MC 210 PetscDTFactorialInt - Compute n! as an integer 211 212 Input Parameter: 213 . n - a non-negative integer 214 215 Output Parameter: 216 . factorial - n! 217 218 Level: beginner 219 220 Note: 221 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. 222 223 .seealso: `PetscDTFactorial()`, `PetscDTBinomialInt()`, `PetscDTBinomial()` 224 M*/ 225 static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial) 226 { 227 PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600}; 228 229 PetscFunctionBegin; 230 *factorial = -1; 231 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); 232 if (n <= 12) { 233 *factorial = facLookup[n]; 234 } else { 235 PetscInt f = facLookup[12]; 236 PetscInt i; 237 238 for (i = 13; i < n + 1; ++i) f *= i; 239 *factorial = f; 240 } 241 PetscFunctionReturn(PETSC_SUCCESS); 242 } 243 244 /*MC 245 PetscDTBinomial - Approximate the binomial coefficient `n` choose `k` 246 247 Input Parameters: 248 + n - a non-negative integer 249 - k - an integer between 0 and `n`, inclusive 250 251 Output Parameter: 252 . binomial - approximation of the binomial coefficient `n` choose `k` 253 254 Level: beginner 255 256 .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()` 257 M*/ 258 static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial) 259 { 260 PetscFunctionBeginHot; 261 *binomial = -1.0; 262 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); 263 if (n <= 3) { 264 PetscInt binomLookup[4][4] = { 265 {1, 0, 0, 0}, 266 {1, 1, 0, 0}, 267 {1, 2, 1, 0}, 268 {1, 3, 3, 1} 269 }; 270 271 *binomial = (PetscReal)binomLookup[n][k]; 272 } else { 273 PetscReal binom = 1.0; 274 275 k = PetscMin(k, n - k); 276 for (PetscInt i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1); 277 *binomial = binom; 278 } 279 PetscFunctionReturn(PETSC_SUCCESS); 280 } 281 282 /*MC 283 PetscDTBinomialInt - Compute the binomial coefficient `n` choose `k` 284 285 Input Parameters: 286 + n - a non-negative integer 287 - k - an integer between 0 and `n`, inclusive 288 289 Output Parameter: 290 . binomial - the binomial coefficient `n` choose `k` 291 292 Level: beginner 293 294 Note: 295 This is limited by integers that can be represented by `PetscInt`. 296 297 Use `PetscDTBinomial()` for real number approximations of larger values 298 299 .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTEnumPerm()` 300 M*/ 301 static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial) 302 { 303 PetscInt bin; 304 305 PetscFunctionBegin; 306 *binomial = -1; 307 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); 308 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); 309 if (n <= 3) { 310 PetscInt binomLookup[4][4] = { 311 {1, 0, 0, 0}, 312 {1, 1, 0, 0}, 313 {1, 2, 1, 0}, 314 {1, 3, 3, 1} 315 }; 316 317 bin = binomLookup[n][k]; 318 } else { 319 PetscInt binom = 1; 320 321 k = PetscMin(k, n - k); 322 for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1); 323 bin = binom; 324 } 325 *binomial = bin; 326 PetscFunctionReturn(PETSC_SUCCESS); 327 } 328 329 /* the following inline routines should be not be inline routines and then Fortran binding can be built automatically */ 330 #define PeOp 331 332 /*MC 333 PetscDTEnumPerm - Get a permutation of `n` integers from its encoding into the integers [0, n!) as a sequence of swaps. 334 335 Input Parameters: 336 + n - a non-negative integer (see note about limits below) 337 - k - an integer in [0, n!) 338 339 Output Parameters: 340 + perm - the permuted list of the integers [0, ..., n-1] 341 - isOdd - if not `NULL`, returns whether the permutation used an even or odd number of swaps. 342 343 Level: intermediate 344 345 Notes: 346 A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation, 347 by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in 348 some position j >= i. This swap is encoded as the difference (j - i). The difference d_i at step i is less than 349 (n - i). This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number 350 (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}. 351 352 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. 353 354 .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTPermIndex()` 355 M*/ 356 static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PeOp PetscBool *isOdd) 357 { 358 PetscInt odd = 0; 359 PetscInt i; 360 PetscInt work[PETSC_FACTORIAL_MAX]; 361 PetscInt *w; 362 363 PetscFunctionBegin; 364 if (isOdd) *isOdd = PETSC_FALSE; 365 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); 366 if (n >= 2) { 367 w = &work[n - 2]; 368 for (i = 2; i <= n; i++) { 369 *(w--) = k % i; 370 k /= i; 371 } 372 } 373 for (i = 0; i < n; i++) perm[i] = i; 374 for (i = 0; i < n - 1; i++) { 375 PetscInt s = work[i]; 376 PetscInt swap = perm[i]; 377 378 perm[i] = perm[i + s]; 379 perm[i + s] = swap; 380 odd ^= (!!s); 381 } 382 if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 383 PetscFunctionReturn(PETSC_SUCCESS); 384 } 385 386 /*MC 387 PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!). This inverts `PetscDTEnumPerm()`. 388 389 Input Parameters: 390 + n - a non-negative integer (see note about limits below) 391 - perm - the permuted list of the integers [0, ..., n-1] 392 393 Output Parameters: 394 + k - an integer in [0, n!) 395 - isOdd - if not `NULL`, returns whether the permutation used an even or odd number of swaps. 396 397 Level: beginner 398 399 Note: 400 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. 401 402 .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()` 403 M*/ 404 static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PeOp PetscBool *isOdd) 405 { 406 PetscInt odd = 0; 407 PetscInt i, idx; 408 PetscInt work[PETSC_FACTORIAL_MAX]; 409 PetscInt iwork[PETSC_FACTORIAL_MAX]; 410 411 PetscFunctionBeginHot; 412 *k = -1; 413 if (isOdd) *isOdd = PETSC_FALSE; 414 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); 415 for (i = 0; i < n; i++) work[i] = i; /* partial permutation */ 416 for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */ 417 for (idx = 0, i = 0; i < n - 1; i++) { 418 PetscInt j = perm[i]; 419 PetscInt icur = work[i]; 420 PetscInt jloc = iwork[j]; 421 PetscInt diff = jloc - i; 422 423 idx = idx * (n - i) + diff; 424 /* swap (i, jloc) */ 425 work[i] = j; 426 work[jloc] = icur; 427 iwork[j] = i; 428 iwork[icur] = jloc; 429 odd ^= (!!diff); 430 } 431 *k = idx; 432 if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 433 PetscFunctionReturn(PETSC_SUCCESS); 434 } 435 436 /*MC 437 PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k). 438 The encoding is in lexicographic order. 439 440 Input Parameters: 441 + n - a non-negative integer (see note about limits below) 442 . k - an integer in [0, n] 443 - j - an index in [0, n choose k) 444 445 Output Parameter: 446 . subset - the jth subset of size k of the integers [0, ..., n - 1] 447 448 Level: beginner 449 450 Note: 451 Limited by arguments such that `n` choose `k` can be represented by `PetscInt` 452 453 .seealso: `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()` 454 M*/ 455 static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset) 456 { 457 PetscInt Nk; 458 459 PetscFunctionBeginHot; 460 PetscCall(PetscDTBinomialInt(n, k, &Nk)); 461 for (PetscInt i = 0, l = 0; i < n && l < k; i++) { 462 PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 463 PetscInt Nminusk = Nk - Nminuskminus; 464 465 if (j < Nminuskminus) { 466 subset[l++] = i; 467 Nk = Nminuskminus; 468 } else { 469 j -= Nminuskminus; 470 Nk = Nminusk; 471 } 472 } 473 PetscFunctionReturn(PETSC_SUCCESS); 474 } 475 476 /*MC 477 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. 478 This is the inverse of `PetscDTEnumSubset`. 479 480 Input Parameters: 481 + n - a non-negative integer (see note about limits below) 482 . k - an integer in [0, n] 483 - subset - an ordered subset of the integers [0, ..., n - 1] 484 485 Output Parameter: 486 . index - the rank of the subset in lexicographic order 487 488 Level: beginner 489 490 Note: 491 Limited by arguments such that `n` choose `k` can be represented by `PetscInt` 492 493 .seealso: `PetscDTEnumSubset()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()` 494 M*/ 495 static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index) 496 { 497 PetscInt j = 0, Nk; 498 499 PetscFunctionBegin; 500 *index = -1; 501 PetscCall(PetscDTBinomialInt(n, k, &Nk)); 502 for (PetscInt i = 0, l = 0; i < n && l < k; i++) { 503 PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 504 PetscInt Nminusk = Nk - Nminuskminus; 505 506 if (subset[l] == i) { 507 l++; 508 Nk = Nminuskminus; 509 } else { 510 j += Nminuskminus; 511 Nk = Nminusk; 512 } 513 } 514 *index = j; 515 PetscFunctionReturn(PETSC_SUCCESS); 516 } 517 518 /*MC 519 PetscDTEnumSplit - 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. 520 521 Input Parameters: 522 + n - a non-negative integer (see note about limits below) 523 . k - an integer in [0, n] 524 - j - an index in [0, n choose k) 525 526 Output Parameters: 527 + perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set. 528 - isOdd - if not `NULL`, return whether perm is an even or odd permutation. 529 530 Level: beginner 531 532 Note: 533 Limited by arguments such that `n` choose `k` can be represented by `PetscInt` 534 535 .seealso: `PetscDTEnumSubset()`, `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, 536 `PetscDTPermIndex()` 537 M*/ 538 static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PeOp PetscBool *isOdd) 539 { 540 PetscInt i, l, m, Nk, odd = 0; 541 PetscInt *subcomp = PetscSafePointerPlusOffset(perm, k); 542 543 PetscFunctionBegin; 544 if (isOdd) *isOdd = PETSC_FALSE; 545 PetscCall(PetscDTBinomialInt(n, k, &Nk)); 546 for (i = 0, l = 0, m = 0; i < n && l < k; i++) { 547 PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 548 PetscInt Nminusk = Nk - Nminuskminus; 549 550 if (j < Nminuskminus) { 551 perm[l++] = i; 552 Nk = Nminuskminus; 553 } else { 554 subcomp[m++] = i; 555 j -= Nminuskminus; 556 odd ^= ((k - l) & 1); 557 Nk = Nminusk; 558 } 559 } 560 for (; i < n; i++) subcomp[m++] = i; 561 if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 562 PetscFunctionReturn(PETSC_SUCCESS); 563 } 564 565 struct _n_PetscTabulation { 566 PetscInt K; /* Indicates a k-jet, namely tabulated derivatives up to order k */ 567 PetscInt Nr; /* The number of tabulation replicas (often 1) */ 568 PetscInt Np; /* The number of tabulation points in a replica */ 569 PetscInt Nb; /* The number of functions tabulated */ 570 PetscInt Nc; /* The number of function components */ 571 PetscInt cdim; /* The coordinate dimension */ 572 PetscReal **T; /* The tabulation T[K] of functions and their derivatives 573 T[0] = B[Nr*Np][Nb][Nc]: The basis function values at quadrature points 574 T[1] = D[Nr*Np][Nb][Nc][cdim]: The basis function derivatives at quadrature points 575 T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */ 576 }; 577 578 /*S 579 PetscTabulation - PETSc object that manages tabulations for finite element methods. 580 581 Level: intermediate 582 583 Note: 584 This is a pointer to a C struct, hence the data in it may be accessed directly. 585 586 Fortran Note: 587 Use `PetscTabulationGetData()` and `PetscTabulationRestoreData()` to access the arrays in the tabulation. 588 589 Developer Note: 590 TODO: put the meaning of the struct fields in this manual page 591 592 .seealso: `PetscTabulationDestroy()`, `PetscFECreateTabulation()`, `PetscFEGetCellTabulation()` 593 S*/ 594 typedef struct _n_PetscTabulation *PetscTabulation; 595 596 typedef PetscErrorCode (*PetscProbFunc)(const PetscReal[], const PetscReal[], PetscReal[]); 597 598 typedef enum { 599 DTPROB_DENSITY_CONSTANT, 600 DTPROB_DENSITY_GAUSSIAN, 601 DTPROB_DENSITY_MAXWELL_BOLTZMANN, 602 DTPROB_NUM_DENSITY 603 } DTProbDensityType; 604 PETSC_EXTERN const char *const DTProbDensityTypes[]; 605 606 PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]); 607 PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]); 608 PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]); 609 PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]); 610 PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]); 611 PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]); 612 PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 613 PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 614 PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 615 PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]); 616 PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]); 617 PETSC_EXTERN PetscErrorCode PetscPDFGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]); 618 PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]); 619 PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 620 PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 621 PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 622 PETSC_EXTERN PetscErrorCode PetscPDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 623 PETSC_EXTERN PetscErrorCode PetscCDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 624 PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 625 PETSC_EXTERN PetscErrorCode PetscPDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 626 PETSC_EXTERN PetscErrorCode PetscCDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 627 PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 628 PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFunc *, PetscProbFunc *, PetscProbFunc *); 629 630 #include <petscvec.h> 631 632 PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFunc, PetscReal *); 633 PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatisticWeighted(Vec, Vec, PetscProbFunc, PetscReal *); 634 PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatisticMagnitude(Vec, PetscProbFunc, PetscReal *); 635