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