xref: /petsc/include/petscdt.h (revision 8cd1e013f4f44ccfc597f967d47b0f55a2d81953)
137045ce4SJed Brown /*
237045ce4SJed Brown   Common tools for constructing discretizations
337045ce4SJed Brown */
426bd1501SBarry Smith #if !defined(PETSCDT_H)
526bd1501SBarry Smith #define PETSCDT_H
637045ce4SJed Brown 
737045ce4SJed Brown #include <petscsys.h>
837045ce4SJed Brown 
921454ff5SMatthew G. Knepley /*S
1021454ff5SMatthew G. Knepley   PetscQuadrature - Quadrature rule for integration.
1121454ff5SMatthew G. Knepley 
12329bbf4eSMatthew G. Knepley   Level: beginner
1321454ff5SMatthew G. Knepley 
1421454ff5SMatthew G. Knepley .seealso:  PetscQuadratureCreate(), PetscQuadratureDestroy()
1521454ff5SMatthew G. Knepley S*/
1621454ff5SMatthew G. Knepley typedef struct _p_PetscQuadrature *PetscQuadrature;
1721454ff5SMatthew G. Knepley 
188272889dSSatish Balay /*E
19916e780bShannah_mairs   PetscGaussLobattoLegendreCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights
208272889dSSatish Balay 
218272889dSSatish Balay   Level: intermediate
228272889dSSatish Balay 
23f2e8fe4dShannah_mairs $  PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA - compute the nodes via linear algebra
24d410ae54Shannah_mairs $  PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON - compute the nodes by solving a nonlinear equation with Newton's method
258272889dSSatish Balay 
268272889dSSatish Balay E*/
27f2e8fe4dShannah_mairs typedef enum {PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON} PetscGaussLobattoLegendreCreateType;
288272889dSSatish Balay 
2921454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *);
30c9638911SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *);
31bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt*);
32bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt);
33a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt*);
34a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt);
35a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt*, PetscInt*, PetscInt*, const PetscReal *[], const PetscReal *[]);
36a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal [], const PetscReal []);
3721454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer);
3821454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *);
39a0845e3aSMatthew G. Knepley 
4089710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *);
4189710940SMatthew G. Knepley 
4237045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*);
4337045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt,PetscReal,PetscReal,PetscReal*,PetscReal*);
44916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt,PetscGaussLobattoLegendreCreateType,PetscReal*,PetscReal*);
45194825f6SJed Brown PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt,PetscInt,const PetscReal*,PetscInt,const PetscReal*,PetscReal*);
46a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
47a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
4837045ce4SJed Brown 
49b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
50b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *);
51d525116cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *);
52b3c0f97bSTom Klotz 
53916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *);
54916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
55916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
56916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
57916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
58916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
59916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
60916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
61916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
62916e780bShannah_mairs 
631a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
641a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
651a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
661a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
671a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
681a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
691a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
70dda711d0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
711a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
721a989b97SToby Isaac 
73fad4db65SToby Isaac #if defined(PETSC_USE_64BIT_INDICES)
74fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 20
75fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX  61
76fad4db65SToby Isaac #else
77fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 12
78fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX  29
79fad4db65SToby Isaac #endif
80fad4db65SToby Isaac 
81fad4db65SToby Isaac /*MC
82fad4db65SToby Isaac    PetscDTFactorial - Approximate n! as a real number
83fad4db65SToby Isaac 
84fad4db65SToby Isaac    Input Arguments:
85fad4db65SToby Isaac 
86fad4db65SToby Isaac .  n - a non-negative integer
87fad4db65SToby Isaac 
88fad4db65SToby Isaac    Output Arguments;
89fad4db65SToby Isaac 
90fad4db65SToby Isaac .  factorial - n!
91fad4db65SToby Isaac 
92fad4db65SToby Isaac    Level: beginner
93fad4db65SToby Isaac M*/
94fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial)
95fad4db65SToby Isaac {
96fad4db65SToby Isaac   PetscReal f = 1.0;
97fad4db65SToby Isaac   PetscInt  i;
98fad4db65SToby Isaac 
99fad4db65SToby Isaac   PetscFunctionBegin;
100fad4db65SToby Isaac   for (i = 1; i < n+1; ++i) f *= i;
101fad4db65SToby Isaac   *factorial = f;
102fad4db65SToby Isaac   PetscFunctionReturn(0);
103fad4db65SToby Isaac }
104fad4db65SToby Isaac 
105fad4db65SToby Isaac /*MC
106fad4db65SToby Isaac    PetscDTFactorialInt - Compute n! as an integer
107fad4db65SToby Isaac 
108fad4db65SToby Isaac    Input Arguments:
109fad4db65SToby Isaac 
110fad4db65SToby Isaac .  n - a non-negative integer
111fad4db65SToby Isaac 
112fad4db65SToby Isaac    Output Arguments;
113fad4db65SToby Isaac 
114fad4db65SToby Isaac .  factorial - n!
115fad4db65SToby Isaac 
116fad4db65SToby Isaac    Level: beginner
117fad4db65SToby Isaac 
118fad4db65SToby Isaac    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.
119fad4db65SToby Isaac M*/
120fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial)
121fad4db65SToby Isaac {
122fad4db65SToby Isaac   PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600};
123fad4db65SToby Isaac 
124fad4db65SToby Isaac   PetscFunctionBeginHot;
125fad4db65SToby Isaac   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);
126fad4db65SToby Isaac   if (n <= 12) {
127fad4db65SToby Isaac     *factorial = facLookup[n];
128fad4db65SToby Isaac   } else {
129fad4db65SToby Isaac     PetscInt f = facLookup[12];
130fad4db65SToby Isaac     PetscInt i;
131fad4db65SToby Isaac 
132fad4db65SToby Isaac     for (i = 13; i < n+1; ++i) f *= i;
133fad4db65SToby Isaac     *factorial = f;
134fad4db65SToby Isaac   }
135fad4db65SToby Isaac   PetscFunctionReturn(0);
136fad4db65SToby Isaac }
137fad4db65SToby Isaac 
138fad4db65SToby Isaac /*MC
139fad4db65SToby Isaac    PetscDTBinomial - Approximate the binomial coefficient "n choose k"
140fad4db65SToby Isaac 
141fad4db65SToby Isaac    Input Arguments:
142fad4db65SToby Isaac 
143fad4db65SToby Isaac +  n - a non-negative integer
144fad4db65SToby Isaac -  k - an integer between 0 and n, inclusive
145fad4db65SToby Isaac 
146fad4db65SToby Isaac    Output Arguments;
147fad4db65SToby Isaac 
148fad4db65SToby Isaac .  binomial - approximation of the binomial coefficient n choose k
149fad4db65SToby Isaac 
150fad4db65SToby Isaac    Level: beginner
151fad4db65SToby Isaac M*/
152fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
1531a989b97SToby Isaac {
1541a989b97SToby Isaac   PetscFunctionBeginHot;
155fad4db65SToby Isaac   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);
1561a989b97SToby Isaac   if (n <= 3) {
1571a989b97SToby Isaac     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};
1581a989b97SToby Isaac 
1591a989b97SToby Isaac     *binomial = binomLookup[n][k];
1601a989b97SToby Isaac   } else {
161fad4db65SToby Isaac     PetscReal binom = 1.;
1621a989b97SToby Isaac     PetscInt  i;
1631a989b97SToby Isaac 
1641a989b97SToby Isaac     k = PetscMin(k, n - k);
1651a989b97SToby Isaac     for (i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
1661a989b97SToby Isaac     *binomial = binom;
1671a989b97SToby Isaac   }
1681a989b97SToby Isaac   PetscFunctionReturn(0);
1691a989b97SToby Isaac }
1701a989b97SToby Isaac 
171fad4db65SToby Isaac /*MC
172fad4db65SToby Isaac    PetscDTBinomialInt - Compute the binomial coefficient "n choose k"
173fad4db65SToby Isaac 
174fad4db65SToby Isaac    Input Arguments:
175fad4db65SToby Isaac 
176fad4db65SToby Isaac +  n - a non-negative integer
177fad4db65SToby Isaac -  k - an integer between 0 and n, inclusive
178fad4db65SToby Isaac 
179fad4db65SToby Isaac    Output Arguments;
180fad4db65SToby Isaac 
181fad4db65SToby Isaac .  binomial - the binomial coefficient n choose k
182fad4db65SToby Isaac 
183fad4db65SToby Isaac    Note: this is limited by integers that can be represented by PetscInt
184fad4db65SToby Isaac 
185fad4db65SToby Isaac    Level: beginner
186fad4db65SToby Isaac M*/
187fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial)
188fad4db65SToby Isaac {
189fad4db65SToby Isaac   PetscFunctionBeginHot;
190fad4db65SToby Isaac   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);
191fad4db65SToby Isaac   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);
192fad4db65SToby Isaac   if (n <= 3) {
193fad4db65SToby Isaac     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};
194fad4db65SToby Isaac 
195fad4db65SToby Isaac     *binomial = binomLookup[n][k];
196fad4db65SToby Isaac   } else {
197fad4db65SToby Isaac     PetscInt  binom = 1;
198fad4db65SToby Isaac     PetscInt  i;
199fad4db65SToby Isaac 
200fad4db65SToby Isaac     k = PetscMin(k, n - k);
201fad4db65SToby Isaac     for (i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
202fad4db65SToby Isaac     *binomial = binom;
203fad4db65SToby Isaac   }
204fad4db65SToby Isaac   PetscFunctionReturn(0);
205fad4db65SToby Isaac }
206fad4db65SToby Isaac 
207fad4db65SToby Isaac /*MC
208fad4db65SToby Isaac    PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps.
209fad4db65SToby Isaac 
210fad4db65SToby Isaac    A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation,
211fad4db65SToby Isaac    by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in
212fad4db65SToby Isaac    some position j >= i.  We encode this swap as the difference (j - i).  The difference d_i at step i is less than
213fad4db65SToby Isaac    (n - i).  We encode this sequence of n-1 differences [d_0, ..., d_{n-2}] as the number
214fad4db65SToby Isaac    (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}.
215fad4db65SToby Isaac 
216fad4db65SToby Isaac    Input Arguments:
217fad4db65SToby Isaac 
218fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
219*8cd1e013SToby Isaac -  k - an integer in [0, n!)
220fad4db65SToby Isaac 
221fad4db65SToby Isaac    Output Arguments:
222fad4db65SToby Isaac 
223fad4db65SToby Isaac +  perm - the permuted list of the integers [0, ..., n-1]
224*8cd1e013SToby Isaac -  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.
225fad4db65SToby Isaac 
226fad4db65SToby Isaac    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.
227fad4db65SToby Isaac 
228fad4db65SToby Isaac    Level: beginner
229fad4db65SToby Isaac M*/
230fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd)
2311a989b97SToby Isaac {
2321a989b97SToby Isaac   PetscInt  odd = 0;
2331a989b97SToby Isaac   PetscInt  i;
234fad4db65SToby Isaac   PetscInt  work[PETSC_FACTORIAL_MAX];
235fad4db65SToby Isaac   PetscInt *w;
2361a989b97SToby Isaac 
2371a989b97SToby Isaac   PetscFunctionBeginHot;
238fad4db65SToby Isaac   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);
239fad4db65SToby Isaac   w = &work[n - 2];
2401a989b97SToby Isaac   for (i = 2; i <= n; i++) {
2411a989b97SToby Isaac     *(w--) = k % i;
2421a989b97SToby Isaac     k /= i;
2431a989b97SToby Isaac   }
2441a989b97SToby Isaac   for (i = 0; i < n; i++) perm[i] = i;
2451a989b97SToby Isaac   for (i = 0; i < n - 1; i++) {
2461a989b97SToby Isaac     PetscInt s = work[i];
2471a989b97SToby Isaac     PetscInt swap = perm[i];
2481a989b97SToby Isaac 
2491a989b97SToby Isaac     perm[i] = perm[i + s];
2501a989b97SToby Isaac     perm[i + s] = swap;
2511a989b97SToby Isaac     odd ^= (!!s);
2521a989b97SToby Isaac   }
2531a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
2541a989b97SToby Isaac   PetscFunctionReturn(0);
2551a989b97SToby Isaac }
2561a989b97SToby Isaac 
257fad4db65SToby Isaac /*MC
258*8cd1e013SToby Isaac    PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!).  This inverts PetscDTEnumPerm.
259*8cd1e013SToby Isaac 
260*8cd1e013SToby Isaac    Input Arguments:
261*8cd1e013SToby Isaac 
262*8cd1e013SToby Isaac +  n - a non-negative integer (see note about limits below)
263*8cd1e013SToby Isaac -  perm - the permuted list of the integers [0, ..., n-1]
264*8cd1e013SToby Isaac 
265*8cd1e013SToby Isaac    Output Arguments:
266*8cd1e013SToby Isaac 
267*8cd1e013SToby Isaac +  k - an integer in [0, n!)
268*8cd1e013SToby Isaac .  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.
269*8cd1e013SToby Isaac 
270*8cd1e013SToby Isaac    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.
271*8cd1e013SToby Isaac 
272*8cd1e013SToby Isaac    Level: beginner
273*8cd1e013SToby Isaac M*/
274*8cd1e013SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd)
275*8cd1e013SToby Isaac {
276*8cd1e013SToby Isaac   PetscInt  odd = 0;
277*8cd1e013SToby Isaac   PetscInt  i, idx;
278*8cd1e013SToby Isaac   PetscInt  work[PETSC_FACTORIAL_MAX];
279*8cd1e013SToby Isaac   PetscInt  iwork[PETSC_FACTORIAL_MAX];
280*8cd1e013SToby Isaac 
281*8cd1e013SToby Isaac   PetscFunctionBeginHot;
282*8cd1e013SToby Isaac   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*8cd1e013SToby Isaac   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
284*8cd1e013SToby Isaac   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
285*8cd1e013SToby Isaac   for (idx = 0, i = 0; i < n - 1; i++) {
286*8cd1e013SToby Isaac     PetscInt j = perm[i];
287*8cd1e013SToby Isaac     PetscInt icur = work[i];
288*8cd1e013SToby Isaac     PetscInt jloc = iwork[j];
289*8cd1e013SToby Isaac     PetscInt diff = jloc - i;
290*8cd1e013SToby Isaac 
291*8cd1e013SToby Isaac     idx = idx * (n - i) + diff;
292*8cd1e013SToby Isaac     /* swap (i, jloc) */
293*8cd1e013SToby Isaac     work[i] = j;
294*8cd1e013SToby Isaac     work[jloc] = icur;
295*8cd1e013SToby Isaac     iwork[j] = i;
296*8cd1e013SToby Isaac     iwork[icur] = jloc;
297*8cd1e013SToby Isaac     odd ^= (!!diff);
298*8cd1e013SToby Isaac   }
299*8cd1e013SToby Isaac   *k = idx;
300*8cd1e013SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
301*8cd1e013SToby Isaac   PetscFunctionReturn(0);
302*8cd1e013SToby Isaac }
303*8cd1e013SToby Isaac 
304*8cd1e013SToby Isaac /*MC
305fad4db65SToby Isaac    PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k).
306fad4db65SToby Isaac    The encoding is in lexicographic order.
307fad4db65SToby Isaac 
308fad4db65SToby Isaac    Input Arguments:
309fad4db65SToby Isaac 
310fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
311fad4db65SToby Isaac .  k - an integer in [0, n]
312fad4db65SToby Isaac -  j - an index in [0, n choose k)
313fad4db65SToby Isaac 
314fad4db65SToby Isaac    Output Arguments:
315fad4db65SToby Isaac 
316fad4db65SToby Isaac .  subset - the jth subset of size k of the integers [0, ..., n - 1]
317fad4db65SToby Isaac 
318fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
319fad4db65SToby Isaac 
320fad4db65SToby Isaac    Level: beginner
321fad4db65SToby Isaac 
322fad4db65SToby Isaac .seealso: PetscDTSubsetIndex()
323fad4db65SToby Isaac M*/
3241a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
3251a989b97SToby Isaac {
3261a989b97SToby Isaac   PetscInt       Nk, i, l;
3271a989b97SToby Isaac   PetscErrorCode ierr;
3281a989b97SToby Isaac 
3291a989b97SToby Isaac   PetscFunctionBeginHot;
330fad4db65SToby Isaac   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
3311a989b97SToby Isaac   for (i = 0, l = 0; i < n && l < k; i++) {
3321a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
3331a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
3341a989b97SToby Isaac 
3351a989b97SToby Isaac     if (j < Nminuskminus) {
3361a989b97SToby Isaac       subset[l++] = i;
3371a989b97SToby Isaac       Nk = Nminuskminus;
3381a989b97SToby Isaac     } else {
3391a989b97SToby Isaac       j -= Nminuskminus;
3401a989b97SToby Isaac       Nk = Nminusk;
3411a989b97SToby Isaac     }
3421a989b97SToby Isaac   }
3431a989b97SToby Isaac   PetscFunctionReturn(0);
3441a989b97SToby Isaac }
3451a989b97SToby Isaac 
346fad4db65SToby Isaac /*MC
347fad4db65SToby Isaac    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.
348fad4db65SToby Isaac 
349fad4db65SToby Isaac    Input Arguments:
350fad4db65SToby Isaac 
351fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
352fad4db65SToby Isaac .  k - an integer in [0, n]
353fad4db65SToby Isaac -  subset - an ordered subset of the integers [0, ..., n - 1]
354fad4db65SToby Isaac 
355fad4db65SToby Isaac    Output Arguments:
356fad4db65SToby Isaac 
357fad4db65SToby Isaac .  index - the rank of the subset in lexicographic order
358fad4db65SToby Isaac 
359fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
360fad4db65SToby Isaac 
361fad4db65SToby Isaac    Level: beginner
362fad4db65SToby Isaac 
363fad4db65SToby Isaac .seealso: PetscDTEnumSubset()
364fad4db65SToby Isaac M*/
3651a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
3661a989b97SToby Isaac {
3671a989b97SToby Isaac   PetscInt       i, j = 0, l, Nk;
3681a989b97SToby Isaac   PetscErrorCode ierr;
3691a989b97SToby Isaac 
3701a989b97SToby Isaac   PetscFunctionBeginHot;
371fad4db65SToby Isaac   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
3721a989b97SToby Isaac   for (i = 0, l = 0; i < n && l < k; i++) {
3731a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
3741a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
3751a989b97SToby Isaac 
3761a989b97SToby Isaac     if (subset[l] == i) {
3771a989b97SToby Isaac       l++;
3781a989b97SToby Isaac       Nk = Nminuskminus;
3791a989b97SToby Isaac     } else {
3801a989b97SToby Isaac       j += Nminuskminus;
3811a989b97SToby Isaac       Nk = Nminusk;
3821a989b97SToby Isaac     }
3831a989b97SToby Isaac   }
3841a989b97SToby Isaac   *index = j;
3851a989b97SToby Isaac   PetscFunctionReturn(0);
3861a989b97SToby Isaac }
3871a989b97SToby Isaac 
3881a989b97SToby Isaac 
389fad4db65SToby Isaac /*MC
390fad4db65SToby Isaac    PetscDTEnumSubset - Split the integers [0, ..., n - 1] into two complementary ordered subsets, the first of size k and beingthe jth in lexicographic order.
391fad4db65SToby Isaac 
392fad4db65SToby Isaac    Input Arguments:
393fad4db65SToby Isaac 
394fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
395fad4db65SToby Isaac .  k - an integer in [0, n]
396fad4db65SToby Isaac -  j - an index in [0, n choose k)
397fad4db65SToby Isaac 
398fad4db65SToby Isaac    Output Arguments:
399fad4db65SToby Isaac 
400fad4db65SToby Isaac +  perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set.
401fad4db65SToby Isaac -  isOdd - if not NULL, return whether the permutation is even or odd.
402fad4db65SToby Isaac 
403fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
404fad4db65SToby Isaac 
405fad4db65SToby Isaac    Level: beginner
406fad4db65SToby Isaac 
407fad4db65SToby Isaac .seealso: PetscDTEnumSubset(), PetscDTSubsetIndex()
408fad4db65SToby Isaac M*/
409fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd)
4101a989b97SToby Isaac {
4111a989b97SToby Isaac   PetscInt       i, l, m, *subcomp, Nk;
4121a989b97SToby Isaac   PetscInt       odd;
4131a989b97SToby Isaac   PetscErrorCode ierr;
4141a989b97SToby Isaac 
4151a989b97SToby Isaac   PetscFunctionBeginHot;
416fad4db65SToby Isaac   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
4171a989b97SToby Isaac   odd = 0;
418fad4db65SToby Isaac   subcomp = &perm[k];
4191a989b97SToby Isaac   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
4201a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
4211a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
4221a989b97SToby Isaac 
4231a989b97SToby Isaac     if (j < Nminuskminus) {
424fad4db65SToby Isaac       perm[l++] = i;
4251a989b97SToby Isaac       Nk = Nminuskminus;
4261a989b97SToby Isaac     } else {
4271a989b97SToby Isaac       subcomp[m++] = i;
4281a989b97SToby Isaac       j -= Nminuskminus;
4291a989b97SToby Isaac       odd ^= ((k - l) & 1);
4301a989b97SToby Isaac       Nk = Nminusk;
4311a989b97SToby Isaac     }
4321a989b97SToby Isaac   }
4331a989b97SToby Isaac   for (; i < n; i++) {
4341a989b97SToby Isaac     subcomp[m++] = i;
4351a989b97SToby Isaac   }
4361a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
4371a989b97SToby Isaac   PetscFunctionReturn(0);
4381a989b97SToby Isaac }
4391a989b97SToby Isaac 
44037045ce4SJed Brown #endif
441