1 /* 2 Common tools for constructing discretizations 3 */ 4 #if !defined(PETSCDT_H) 5 #define PETSCDT_H 6 7 #include <petscsys.h> 8 9 /*S 10 PetscQuadrature - Quadrature rule for integration. 11 12 Level: beginner 13 14 .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy() 15 S*/ 16 typedef struct _p_PetscQuadrature *PetscQuadrature; 17 18 /*E 19 PetscGaussLobattoLegendreCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights 20 21 Level: intermediate 22 23 $ PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA - compute the nodes via linear algebra 24 $ PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON - compute the nodes by solving a nonlinear equation with Newton's method 25 26 E*/ 27 typedef enum {PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON} PetscGaussLobattoLegendreCreateType; 28 29 PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *); 30 PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *); 31 PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt*); 32 PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt); 33 PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt*); 34 PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt); 35 PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt*, PetscInt*, PetscInt*, const PetscReal *[], const PetscReal *[]); 36 PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal [], const PetscReal []); 37 PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer); 38 PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *); 39 40 PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *); 41 42 PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *); 43 44 PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*); 45 PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt,PetscReal,PetscReal,PetscReal*,PetscReal*); 46 PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt,PetscGaussLobattoLegendreCreateType,PetscReal*,PetscReal*); 47 PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt,PetscInt,const PetscReal*,PetscInt,const PetscReal*,PetscReal*); 48 PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*); 49 PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*); 50 51 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *); 52 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *); 53 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *); 54 55 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *); 56 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 57 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 58 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 59 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 60 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 61 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 62 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 63 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 64 65 PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 66 PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 67 PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 68 PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *); 69 PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *); 70 PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 71 PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *); 72 PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]); 73 PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 74 75 #if defined(PETSC_USE_64BIT_INDICES) 76 #define PETSC_FACTORIAL_MAX 20 77 #define PETSC_BINOMIAL_MAX 61 78 #else 79 #define PETSC_FACTORIAL_MAX 12 80 #define PETSC_BINOMIAL_MAX 29 81 #endif 82 83 /*MC 84 PetscDTFactorial - Approximate n! as a real number 85 86 Input Arguments: 87 . n - a non-negative integer 88 89 Output Arguments: 90 . factorial - n! 91 92 Level: beginner 93 M*/ 94 PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial) 95 { 96 PetscReal f = 1.0; 97 PetscInt i; 98 99 PetscFunctionBegin; 100 *factorial = -1.; 101 if (n < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %D\n", n); 102 for (i = 1; i < n+1; ++i) f *= i; 103 *factorial = f; 104 PetscFunctionReturn(0); 105 } 106 107 /*MC 108 PetscDTFactorialInt - Compute n! as an integer 109 110 Input Arguments: 111 . n - a non-negative integer 112 113 Output Arguments: 114 . factorial - n! 115 116 Level: beginner 117 118 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. 119 M*/ 120 PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial) 121 { 122 PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600}; 123 124 PetscFunctionBegin; 125 *factorial = -1; 126 if (n < 0 || n > PETSC_FACTORIAL_MAX) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of elements %D is not in supported range [0,%D]\n",n,PETSC_FACTORIAL_MAX); 127 if (n <= 12) { 128 *factorial = facLookup[n]; 129 } else { 130 PetscInt f = facLookup[12]; 131 PetscInt i; 132 133 for (i = 13; i < n+1; ++i) f *= i; 134 *factorial = f; 135 } 136 PetscFunctionReturn(0); 137 } 138 139 /*MC 140 PetscDTBinomial - Approximate the binomial coefficient "n choose k" 141 142 Input Arguments: 143 + n - a non-negative integer 144 - k - an integer between 0 and n, inclusive 145 146 Output Arguments: 147 . binomial - approximation of the binomial coefficient n choose k 148 149 Level: beginner 150 M*/ 151 PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial) 152 { 153 PetscFunctionBeginHot; 154 if (n < 0 || k < 0 || k > n) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial arguments (%D %D) must be non-negative, k <= n\n", n, k); 155 if (n <= 3) { 156 PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}}; 157 158 *binomial = binomLookup[n][k]; 159 } else { 160 PetscReal binom = 1.; 161 PetscInt i; 162 163 k = PetscMin(k, n - k); 164 for (i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1); 165 *binomial = binom; 166 } 167 PetscFunctionReturn(0); 168 } 169 170 /*MC 171 PetscDTBinomialInt - Compute the binomial coefficient "n choose k" 172 173 Input Arguments: 174 + n - a non-negative integer 175 - k - an integer between 0 and n, inclusive 176 177 Output Arguments: 178 . binomial - the binomial coefficient n choose k 179 180 Note: this is limited by integers that can be represented by PetscInt 181 182 Level: beginner 183 M*/ 184 PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial) 185 { 186 PetscInt bin; 187 188 PetscFunctionBegin; 189 *binomial = -1; 190 if (n < 0 || k < 0 || k > n) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial arguments (%D %D) must be non-negative, k <= n\n", n, k); 191 if (n > PETSC_BINOMIAL_MAX) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial elements %D is larger than max for PetscInt, %D\n", n, PETSC_BINOMIAL_MAX); 192 if (n <= 3) { 193 PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}}; 194 195 bin = binomLookup[n][k]; 196 } else { 197 PetscInt binom = 1; 198 PetscInt i; 199 200 k = PetscMin(k, n - k); 201 for (i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1); 202 bin = binom; 203 } 204 *binomial = bin; 205 PetscFunctionReturn(0); 206 } 207 208 /*MC 209 PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps. 210 211 A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation, 212 by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in 213 some position j >= i. This swap is encoded as the difference (j - i). The difference d_i at step i is less than 214 (n - i). This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number 215 (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}. 216 217 Input Arguments: 218 + n - a non-negative integer (see note about limits below) 219 - k - an integer in [0, n!) 220 221 Output Arguments: 222 + perm - the permuted list of the integers [0, ..., n-1] 223 - isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps. 224 225 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. 226 227 Level: beginner 228 M*/ 229 PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd) 230 { 231 PetscInt odd = 0; 232 PetscInt i; 233 PetscInt work[PETSC_FACTORIAL_MAX]; 234 PetscInt *w; 235 236 PetscFunctionBegin; 237 if (isOdd) *isOdd = PETSC_FALSE; 238 if (n < 0 || n > PETSC_FACTORIAL_MAX) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of elements %D is not in supported range [0,%D]\n",n,PETSC_FACTORIAL_MAX); 239 w = &work[n - 2]; 240 for (i = 2; i <= n; i++) { 241 *(w--) = k % i; 242 k /= i; 243 } 244 for (i = 0; i < n; i++) perm[i] = i; 245 for (i = 0; i < n - 1; i++) { 246 PetscInt s = work[i]; 247 PetscInt swap = perm[i]; 248 249 perm[i] = perm[i + s]; 250 perm[i + s] = swap; 251 odd ^= (!!s); 252 } 253 if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 254 PetscFunctionReturn(0); 255 } 256 257 /*MC 258 PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!). This inverts PetscDTEnumPerm. 259 260 Input Arguments: 261 + n - a non-negative integer (see note about limits below) 262 - perm - the permuted list of the integers [0, ..., n-1] 263 264 Output Arguments: 265 + k - an integer in [0, n!) 266 . isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps. 267 268 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. 269 270 Level: beginner 271 M*/ 272 PETSC_STATIC_INLINE PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd) 273 { 274 PetscInt odd = 0; 275 PetscInt i, idx; 276 PetscInt work[PETSC_FACTORIAL_MAX]; 277 PetscInt iwork[PETSC_FACTORIAL_MAX]; 278 279 PetscFunctionBeginHot; 280 *k = -1; 281 if (isOdd) *isOdd = PETSC_FALSE; 282 if (n < 0 || n > PETSC_FACTORIAL_MAX) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of elements %D is not in supported range [0,%D]\n",n,PETSC_FACTORIAL_MAX); 283 for (i = 0; i < n; i++) work[i] = i; /* partial permutation */ 284 for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */ 285 for (idx = 0, i = 0; i < n - 1; i++) { 286 PetscInt j = perm[i]; 287 PetscInt icur = work[i]; 288 PetscInt jloc = iwork[j]; 289 PetscInt diff = jloc - i; 290 291 idx = idx * (n - i) + diff; 292 /* swap (i, jloc) */ 293 work[i] = j; 294 work[jloc] = icur; 295 iwork[j] = i; 296 iwork[icur] = jloc; 297 odd ^= (!!diff); 298 } 299 *k = idx; 300 if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 301 PetscFunctionReturn(0); 302 } 303 304 /*MC 305 PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k). 306 The encoding is in lexicographic order. 307 308 Input Arguments: 309 + n - a non-negative integer (see note about limits below) 310 . k - an integer in [0, n] 311 - j - an index in [0, n choose k) 312 313 Output Arguments: 314 . subset - the jth subset of size k of the integers [0, ..., n - 1] 315 316 Note: this is limited by arguments such that n choose k can be represented by PetscInt 317 318 Level: beginner 319 320 .seealso: PetscDTSubsetIndex() 321 M*/ 322 PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset) 323 { 324 PetscInt Nk, i, l; 325 PetscErrorCode ierr; 326 327 PetscFunctionBeginHot; 328 ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr); 329 for (i = 0, l = 0; i < n && l < k; i++) { 330 PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 331 PetscInt Nminusk = Nk - Nminuskminus; 332 333 if (j < Nminuskminus) { 334 subset[l++] = i; 335 Nk = Nminuskminus; 336 } else { 337 j -= Nminuskminus; 338 Nk = Nminusk; 339 } 340 } 341 PetscFunctionReturn(0); 342 } 343 344 /*MC 345 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. 346 347 Input Arguments: 348 + n - a non-negative integer (see note about limits below) 349 . k - an integer in [0, n] 350 - subset - an ordered subset of the integers [0, ..., n - 1] 351 352 Output Arguments: 353 . index - the rank of the subset in lexicographic order 354 355 Note: this is limited by arguments such that n choose k can be represented by PetscInt 356 357 Level: beginner 358 359 .seealso: PetscDTEnumSubset() 360 M*/ 361 PETSC_STATIC_INLINE PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index) 362 { 363 PetscInt i, j = 0, l, Nk; 364 PetscErrorCode ierr; 365 366 PetscFunctionBegin; 367 *index = -1; 368 ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr); 369 for (i = 0, l = 0; i < n && l < k; i++) { 370 PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 371 PetscInt Nminusk = Nk - Nminuskminus; 372 373 if (subset[l] == i) { 374 l++; 375 Nk = Nminuskminus; 376 } else { 377 j += Nminuskminus; 378 Nk = Nminusk; 379 } 380 } 381 *index = j; 382 PetscFunctionReturn(0); 383 } 384 385 /*MC 386 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. 387 388 Input Arguments: 389 + n - a non-negative integer (see note about limits below) 390 . k - an integer in [0, n] 391 - j - an index in [0, n choose k) 392 393 Output Arguments: 394 + perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set. 395 - isOdd - if not NULL, return whether perm is an even or odd permutation. 396 397 Note: this is limited by arguments such that n choose k can be represented by PetscInt 398 399 Level: beginner 400 401 .seealso: PetscDTEnumSubset(), PetscDTSubsetIndex() 402 M*/ 403 PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd) 404 { 405 PetscInt i, l, m, *subcomp, Nk; 406 PetscInt odd; 407 PetscErrorCode ierr; 408 409 PetscFunctionBegin; 410 if (isOdd) *isOdd = PETSC_FALSE; 411 ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr); 412 odd = 0; 413 subcomp = &perm[k]; 414 for (i = 0, l = 0, m = 0; i < n && l < k; i++) { 415 PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 416 PetscInt Nminusk = Nk - Nminuskminus; 417 418 if (j < Nminuskminus) { 419 perm[l++] = i; 420 Nk = Nminuskminus; 421 } else { 422 subcomp[m++] = i; 423 j -= Nminuskminus; 424 odd ^= ((k - l) & 1); 425 Nk = Nminusk; 426 } 427 } 428 for (; i < n; i++) { 429 subcomp[m++] = i; 430 } 431 if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 432 PetscFunctionReturn(0); 433 } 434 435 #endif 436