xref: /petsc/include/petscdt.h (revision 28222859993257eecc81314090ee1fbecec5d4dc)
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 
42907761f8SToby Isaac PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *);
43907761f8SToby Isaac 
4437045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*);
4537045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt,PetscReal,PetscReal,PetscReal*,PetscReal*);
46916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt,PetscGaussLobattoLegendreCreateType,PetscReal*,PetscReal*);
47194825f6SJed Brown PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt,PetscInt,const PetscReal*,PetscInt,const PetscReal*,PetscReal*);
48a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
49a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
5037045ce4SJed Brown 
51b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
52b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *);
53d525116cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *);
54b3c0f97bSTom Klotz 
55916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *);
56916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
57916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
58916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
59916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
60916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
61916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
62916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
63916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
64916e780bShannah_mairs 
651a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
661a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
671a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
681a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
691a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
701a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
711a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
72dda711d0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
731a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
741a989b97SToby Isaac 
75fad4db65SToby Isaac #if defined(PETSC_USE_64BIT_INDICES)
76fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 20
77fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX  61
78fad4db65SToby Isaac #else
79fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 12
80fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX  29
81fad4db65SToby Isaac #endif
82fad4db65SToby Isaac 
83fad4db65SToby Isaac /*MC
84fad4db65SToby Isaac    PetscDTFactorial - Approximate n! as a real number
85fad4db65SToby Isaac 
86fad4db65SToby Isaac    Input Arguments:
87fad4db65SToby Isaac .  n - a non-negative integer
88fad4db65SToby Isaac 
89*28222859SToby Isaac    Output Arguments:
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;
100*28222859SToby Isaac   *factorial = -1.;
101*28222859SToby Isaac   if (n < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %D\n", n);
102fad4db65SToby Isaac   for (i = 1; i < n+1; ++i) f *= i;
103fad4db65SToby Isaac   *factorial = f;
104fad4db65SToby Isaac   PetscFunctionReturn(0);
105fad4db65SToby Isaac }
106fad4db65SToby Isaac 
107fad4db65SToby Isaac /*MC
108fad4db65SToby Isaac    PetscDTFactorialInt - Compute n! as an integer
109fad4db65SToby Isaac 
110fad4db65SToby Isaac    Input Arguments:
111fad4db65SToby Isaac .  n - a non-negative integer
112fad4db65SToby Isaac 
113*28222859SToby Isaac    Output Arguments:
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 
124*28222859SToby Isaac   PetscFunctionBegin;
125*28222859SToby Isaac   *factorial = -1;
126fad4db65SToby 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);
127fad4db65SToby Isaac   if (n <= 12) {
128fad4db65SToby Isaac     *factorial = facLookup[n];
129fad4db65SToby Isaac   } else {
130fad4db65SToby Isaac     PetscInt f = facLookup[12];
131fad4db65SToby Isaac     PetscInt i;
132fad4db65SToby Isaac 
133fad4db65SToby Isaac     for (i = 13; i < n+1; ++i) f *= i;
134fad4db65SToby Isaac     *factorial = f;
135fad4db65SToby Isaac   }
136fad4db65SToby Isaac   PetscFunctionReturn(0);
137fad4db65SToby Isaac }
138fad4db65SToby Isaac 
139fad4db65SToby Isaac /*MC
140fad4db65SToby Isaac    PetscDTBinomial - Approximate the binomial coefficient "n choose k"
141fad4db65SToby Isaac 
142fad4db65SToby Isaac    Input Arguments:
143fad4db65SToby Isaac +  n - a non-negative integer
144fad4db65SToby Isaac -  k - an integer between 0 and n, inclusive
145fad4db65SToby Isaac 
146*28222859SToby Isaac    Output Arguments:
147fad4db65SToby Isaac .  binomial - approximation of the binomial coefficient n choose k
148fad4db65SToby Isaac 
149fad4db65SToby Isaac    Level: beginner
150fad4db65SToby Isaac M*/
151fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
1521a989b97SToby Isaac {
1531a989b97SToby Isaac   PetscFunctionBeginHot;
154fad4db65SToby 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);
1551a989b97SToby Isaac   if (n <= 3) {
1561a989b97SToby Isaac     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};
1571a989b97SToby Isaac 
1581a989b97SToby Isaac     *binomial = binomLookup[n][k];
1591a989b97SToby Isaac   } else {
160fad4db65SToby Isaac     PetscReal binom = 1.;
1611a989b97SToby Isaac     PetscInt  i;
1621a989b97SToby Isaac 
1631a989b97SToby Isaac     k = PetscMin(k, n - k);
1641a989b97SToby Isaac     for (i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
1651a989b97SToby Isaac     *binomial = binom;
1661a989b97SToby Isaac   }
1671a989b97SToby Isaac   PetscFunctionReturn(0);
1681a989b97SToby Isaac }
1691a989b97SToby Isaac 
170fad4db65SToby Isaac /*MC
171fad4db65SToby Isaac    PetscDTBinomialInt - Compute the binomial coefficient "n choose k"
172fad4db65SToby Isaac 
173fad4db65SToby Isaac    Input Arguments:
174fad4db65SToby Isaac +  n - a non-negative integer
175fad4db65SToby Isaac -  k - an integer between 0 and n, inclusive
176fad4db65SToby Isaac 
177*28222859SToby Isaac    Output Arguments:
178fad4db65SToby Isaac .  binomial - the binomial coefficient n choose k
179fad4db65SToby Isaac 
180fad4db65SToby Isaac    Note: this is limited by integers that can be represented by PetscInt
181fad4db65SToby Isaac 
182fad4db65SToby Isaac    Level: beginner
183fad4db65SToby Isaac M*/
184fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial)
185fad4db65SToby Isaac {
186*28222859SToby Isaac   PetscInt bin;
187*28222859SToby Isaac 
188*28222859SToby Isaac   PetscFunctionBegin;
189*28222859SToby Isaac   *binomial = -1;
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 
195*28222859SToby Isaac     bin = 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);
202*28222859SToby Isaac     bin = binom;
203fad4db65SToby Isaac   }
204*28222859SToby Isaac   *binomial = bin;
205fad4db65SToby Isaac   PetscFunctionReturn(0);
206fad4db65SToby Isaac }
207fad4db65SToby Isaac 
208fad4db65SToby Isaac /*MC
209fad4db65SToby Isaac    PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps.
210fad4db65SToby Isaac 
211fad4db65SToby Isaac    A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation,
212fad4db65SToby Isaac    by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in
213*28222859SToby Isaac    some position j >= i.  This swap is encoded as the difference (j - i).  The difference d_i at step i is less than
214*28222859SToby Isaac    (n - i).  This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number
215fad4db65SToby Isaac    (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}.
216fad4db65SToby Isaac 
217fad4db65SToby Isaac    Input Arguments:
218fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
2198cd1e013SToby Isaac -  k - an integer in [0, n!)
220fad4db65SToby Isaac 
221fad4db65SToby Isaac    Output Arguments:
222fad4db65SToby Isaac +  perm - the permuted list of the integers [0, ..., n-1]
2238cd1e013SToby Isaac -  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.
224fad4db65SToby Isaac 
225fad4db65SToby 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.
226fad4db65SToby Isaac 
227fad4db65SToby Isaac    Level: beginner
228fad4db65SToby Isaac M*/
229fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd)
2301a989b97SToby Isaac {
2311a989b97SToby Isaac   PetscInt  odd = 0;
2321a989b97SToby Isaac   PetscInt  i;
233fad4db65SToby Isaac   PetscInt  work[PETSC_FACTORIAL_MAX];
234fad4db65SToby Isaac   PetscInt *w;
2351a989b97SToby Isaac 
236*28222859SToby Isaac   PetscFunctionBegin;
237*28222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
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
2588cd1e013SToby Isaac    PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!).  This inverts PetscDTEnumPerm.
2598cd1e013SToby Isaac 
2608cd1e013SToby Isaac    Input Arguments:
2618cd1e013SToby Isaac +  n - a non-negative integer (see note about limits below)
2628cd1e013SToby Isaac -  perm - the permuted list of the integers [0, ..., n-1]
2638cd1e013SToby Isaac 
2648cd1e013SToby Isaac    Output Arguments:
2658cd1e013SToby Isaac +  k - an integer in [0, n!)
2668cd1e013SToby Isaac .  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.
2678cd1e013SToby Isaac 
2688cd1e013SToby 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.
2698cd1e013SToby Isaac 
2708cd1e013SToby Isaac    Level: beginner
2718cd1e013SToby Isaac M*/
2728cd1e013SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd)
2738cd1e013SToby Isaac {
2748cd1e013SToby Isaac   PetscInt  odd = 0;
2758cd1e013SToby Isaac   PetscInt  i, idx;
2768cd1e013SToby Isaac   PetscInt  work[PETSC_FACTORIAL_MAX];
2778cd1e013SToby Isaac   PetscInt  iwork[PETSC_FACTORIAL_MAX];
2788cd1e013SToby Isaac 
2798cd1e013SToby Isaac   PetscFunctionBeginHot;
280*28222859SToby Isaac   *k = -1;
281*28222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
2828cd1e013SToby 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);
2838cd1e013SToby Isaac   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
2848cd1e013SToby Isaac   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
2858cd1e013SToby Isaac   for (idx = 0, i = 0; i < n - 1; i++) {
2868cd1e013SToby Isaac     PetscInt j = perm[i];
2878cd1e013SToby Isaac     PetscInt icur = work[i];
2888cd1e013SToby Isaac     PetscInt jloc = iwork[j];
2898cd1e013SToby Isaac     PetscInt diff = jloc - i;
2908cd1e013SToby Isaac 
2918cd1e013SToby Isaac     idx = idx * (n - i) + diff;
2928cd1e013SToby Isaac     /* swap (i, jloc) */
2938cd1e013SToby Isaac     work[i] = j;
2948cd1e013SToby Isaac     work[jloc] = icur;
2958cd1e013SToby Isaac     iwork[j] = i;
2968cd1e013SToby Isaac     iwork[icur] = jloc;
2978cd1e013SToby Isaac     odd ^= (!!diff);
2988cd1e013SToby Isaac   }
2998cd1e013SToby Isaac   *k = idx;
3008cd1e013SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
3018cd1e013SToby Isaac   PetscFunctionReturn(0);
3028cd1e013SToby Isaac }
3038cd1e013SToby Isaac 
3048cd1e013SToby 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 +  n - a non-negative integer (see note about limits below)
310fad4db65SToby Isaac .  k - an integer in [0, n]
311fad4db65SToby Isaac -  j - an index in [0, n choose k)
312fad4db65SToby Isaac 
313fad4db65SToby Isaac    Output Arguments:
314fad4db65SToby Isaac .  subset - the jth subset of size k of the integers [0, ..., n - 1]
315fad4db65SToby Isaac 
316fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
317fad4db65SToby Isaac 
318fad4db65SToby Isaac    Level: beginner
319fad4db65SToby Isaac 
320fad4db65SToby Isaac .seealso: PetscDTSubsetIndex()
321fad4db65SToby Isaac M*/
3221a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
3231a989b97SToby Isaac {
3241a989b97SToby Isaac   PetscInt       Nk, i, l;
3251a989b97SToby Isaac   PetscErrorCode ierr;
3261a989b97SToby Isaac 
3271a989b97SToby Isaac   PetscFunctionBeginHot;
328fad4db65SToby Isaac   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
3291a989b97SToby Isaac   for (i = 0, l = 0; i < n && l < k; i++) {
3301a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
3311a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
3321a989b97SToby Isaac 
3331a989b97SToby Isaac     if (j < Nminuskminus) {
3341a989b97SToby Isaac       subset[l++] = i;
3351a989b97SToby Isaac       Nk = Nminuskminus;
3361a989b97SToby Isaac     } else {
3371a989b97SToby Isaac       j -= Nminuskminus;
3381a989b97SToby Isaac       Nk = Nminusk;
3391a989b97SToby Isaac     }
3401a989b97SToby Isaac   }
3411a989b97SToby Isaac   PetscFunctionReturn(0);
3421a989b97SToby Isaac }
3431a989b97SToby Isaac 
344fad4db65SToby Isaac /*MC
345fad4db65SToby 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.
346fad4db65SToby Isaac 
347fad4db65SToby Isaac    Input Arguments:
348fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
349fad4db65SToby Isaac .  k - an integer in [0, n]
350fad4db65SToby Isaac -  subset - an ordered subset of the integers [0, ..., n - 1]
351fad4db65SToby Isaac 
352fad4db65SToby Isaac    Output Arguments:
353fad4db65SToby Isaac .  index - the rank of the subset in lexicographic order
354fad4db65SToby Isaac 
355fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
356fad4db65SToby Isaac 
357fad4db65SToby Isaac    Level: beginner
358fad4db65SToby Isaac 
359fad4db65SToby Isaac .seealso: PetscDTEnumSubset()
360fad4db65SToby Isaac M*/
3611a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
3621a989b97SToby Isaac {
3631a989b97SToby Isaac   PetscInt       i, j = 0, l, Nk;
3641a989b97SToby Isaac   PetscErrorCode ierr;
3651a989b97SToby Isaac 
366*28222859SToby Isaac   PetscFunctionBegin;
367*28222859SToby Isaac   *index = -1;
368fad4db65SToby Isaac   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
3691a989b97SToby Isaac   for (i = 0, l = 0; i < n && l < k; i++) {
3701a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
3711a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
3721a989b97SToby Isaac 
3731a989b97SToby Isaac     if (subset[l] == i) {
3741a989b97SToby Isaac       l++;
3751a989b97SToby Isaac       Nk = Nminuskminus;
3761a989b97SToby Isaac     } else {
3771a989b97SToby Isaac       j += Nminuskminus;
3781a989b97SToby Isaac       Nk = Nminusk;
3791a989b97SToby Isaac     }
3801a989b97SToby Isaac   }
3811a989b97SToby Isaac   *index = j;
3821a989b97SToby Isaac   PetscFunctionReturn(0);
3831a989b97SToby Isaac }
3841a989b97SToby Isaac 
385fad4db65SToby Isaac /*MC
386*28222859SToby Isaac    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.
387fad4db65SToby Isaac 
388fad4db65SToby Isaac    Input Arguments:
389fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
390fad4db65SToby Isaac .  k - an integer in [0, n]
391fad4db65SToby Isaac -  j - an index in [0, n choose k)
392fad4db65SToby Isaac 
393fad4db65SToby Isaac    Output Arguments:
394fad4db65SToby Isaac +  perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set.
395*28222859SToby Isaac -  isOdd - if not NULL, return whether perm is an even or odd permutation.
396fad4db65SToby Isaac 
397fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
398fad4db65SToby Isaac 
399fad4db65SToby Isaac    Level: beginner
400fad4db65SToby Isaac 
401fad4db65SToby Isaac .seealso: PetscDTEnumSubset(), PetscDTSubsetIndex()
402fad4db65SToby Isaac M*/
403fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd)
4041a989b97SToby Isaac {
4051a989b97SToby Isaac   PetscInt       i, l, m, *subcomp, Nk;
4061a989b97SToby Isaac   PetscInt       odd;
4071a989b97SToby Isaac   PetscErrorCode ierr;
4081a989b97SToby Isaac 
409*28222859SToby Isaac   PetscFunctionBegin;
410*28222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
411fad4db65SToby Isaac   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
4121a989b97SToby Isaac   odd = 0;
413fad4db65SToby Isaac   subcomp = &perm[k];
4141a989b97SToby Isaac   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
4151a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
4161a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
4171a989b97SToby Isaac 
4181a989b97SToby Isaac     if (j < Nminuskminus) {
419fad4db65SToby Isaac       perm[l++] = i;
4201a989b97SToby Isaac       Nk = Nminuskminus;
4211a989b97SToby Isaac     } else {
4221a989b97SToby Isaac       subcomp[m++] = i;
4231a989b97SToby Isaac       j -= Nminuskminus;
4241a989b97SToby Isaac       odd ^= ((k - l) & 1);
4251a989b97SToby Isaac       Nk = Nminusk;
4261a989b97SToby Isaac     }
4271a989b97SToby Isaac   }
4281a989b97SToby Isaac   for (; i < n; i++) {
4291a989b97SToby Isaac     subcomp[m++] = i;
4301a989b97SToby Isaac   }
4311a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
4321a989b97SToby Isaac   PetscFunctionReturn(0);
4331a989b97SToby Isaac }
4341a989b97SToby Isaac 
43537045ce4SJed Brown #endif
436