xref: /petsc/include/petscdt.h (revision fad4db65a981d86daeaa3759591e8b6204fd5392)
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 
73*fad4db65SToby Isaac #if defined(PETSC_USE_64BIT_INDICES)
74*fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 20
75*fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX  61
76*fad4db65SToby Isaac #else
77*fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 12
78*fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX  29
79*fad4db65SToby Isaac #endif
80*fad4db65SToby Isaac 
81*fad4db65SToby Isaac /*MC
82*fad4db65SToby Isaac    PetscDTFactorial - Approximate n! as a real number
83*fad4db65SToby Isaac 
84*fad4db65SToby Isaac    Input Arguments:
85*fad4db65SToby Isaac 
86*fad4db65SToby Isaac .  n - a non-negative integer
87*fad4db65SToby Isaac 
88*fad4db65SToby Isaac    Output Arguments;
89*fad4db65SToby Isaac 
90*fad4db65SToby Isaac .  factorial - n!
91*fad4db65SToby Isaac 
92*fad4db65SToby Isaac    Level: beginner
93*fad4db65SToby Isaac M*/
94*fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial)
95*fad4db65SToby Isaac {
96*fad4db65SToby Isaac   PetscReal f = 1.0;
97*fad4db65SToby Isaac   PetscInt  i;
98*fad4db65SToby Isaac 
99*fad4db65SToby Isaac   PetscFunctionBegin;
100*fad4db65SToby Isaac   for (i = 1; i < n+1; ++i) f *= i;
101*fad4db65SToby Isaac   *factorial = f;
102*fad4db65SToby Isaac   PetscFunctionReturn(0);
103*fad4db65SToby Isaac }
104*fad4db65SToby Isaac 
105*fad4db65SToby Isaac /*MC
106*fad4db65SToby Isaac    PetscDTFactorialInt - Compute n! as an integer
107*fad4db65SToby Isaac 
108*fad4db65SToby Isaac    Input Arguments:
109*fad4db65SToby Isaac 
110*fad4db65SToby Isaac .  n - a non-negative integer
111*fad4db65SToby Isaac 
112*fad4db65SToby Isaac    Output Arguments;
113*fad4db65SToby Isaac 
114*fad4db65SToby Isaac .  factorial - n!
115*fad4db65SToby Isaac 
116*fad4db65SToby Isaac    Level: beginner
117*fad4db65SToby Isaac 
118*fad4db65SToby 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.
119*fad4db65SToby Isaac M*/
120*fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial)
121*fad4db65SToby Isaac {
122*fad4db65SToby Isaac   PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600};
123*fad4db65SToby Isaac 
124*fad4db65SToby Isaac   PetscFunctionBeginHot;
125*fad4db65SToby 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);
126*fad4db65SToby Isaac   if (n <= 12) {
127*fad4db65SToby Isaac     *factorial = facLookup[n];
128*fad4db65SToby Isaac   } else {
129*fad4db65SToby Isaac     PetscInt f = facLookup[12];
130*fad4db65SToby Isaac     PetscInt i;
131*fad4db65SToby Isaac 
132*fad4db65SToby Isaac     for (i = 13; i < n+1; ++i) f *= i;
133*fad4db65SToby Isaac     *factorial = f;
134*fad4db65SToby Isaac   }
135*fad4db65SToby Isaac   PetscFunctionReturn(0);
136*fad4db65SToby Isaac }
137*fad4db65SToby Isaac 
138*fad4db65SToby Isaac /*MC
139*fad4db65SToby Isaac    PetscDTBinomial - Approximate the binomial coefficient "n choose k"
140*fad4db65SToby Isaac 
141*fad4db65SToby Isaac    Input Arguments:
142*fad4db65SToby Isaac 
143*fad4db65SToby Isaac +  n - a non-negative integer
144*fad4db65SToby Isaac -  k - an integer between 0 and n, inclusive
145*fad4db65SToby Isaac 
146*fad4db65SToby Isaac    Output Arguments;
147*fad4db65SToby Isaac 
148*fad4db65SToby Isaac .  binomial - approximation of the binomial coefficient n choose k
149*fad4db65SToby Isaac 
150*fad4db65SToby Isaac    Level: beginner
151*fad4db65SToby Isaac M*/
152*fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
1531a989b97SToby Isaac {
1541a989b97SToby Isaac   PetscFunctionBeginHot;
155*fad4db65SToby 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 {
161*fad4db65SToby 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 
171*fad4db65SToby Isaac /*MC
172*fad4db65SToby Isaac    PetscDTBinomialInt - Compute the binomial coefficient "n choose k"
173*fad4db65SToby Isaac 
174*fad4db65SToby Isaac    Input Arguments:
175*fad4db65SToby Isaac 
176*fad4db65SToby Isaac +  n - a non-negative integer
177*fad4db65SToby Isaac -  k - an integer between 0 and n, inclusive
178*fad4db65SToby Isaac 
179*fad4db65SToby Isaac    Output Arguments;
180*fad4db65SToby Isaac 
181*fad4db65SToby Isaac .  binomial - the binomial coefficient n choose k
182*fad4db65SToby Isaac 
183*fad4db65SToby Isaac    Note: this is limited by integers that can be represented by PetscInt
184*fad4db65SToby Isaac 
185*fad4db65SToby Isaac    Level: beginner
186*fad4db65SToby Isaac M*/
187*fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial)
188*fad4db65SToby Isaac {
189*fad4db65SToby Isaac   PetscFunctionBeginHot;
190*fad4db65SToby 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);
191*fad4db65SToby 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);
192*fad4db65SToby Isaac   if (n <= 3) {
193*fad4db65SToby Isaac     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};
194*fad4db65SToby Isaac 
195*fad4db65SToby Isaac     *binomial = binomLookup[n][k];
196*fad4db65SToby Isaac   } else {
197*fad4db65SToby Isaac     PetscInt  binom = 1;
198*fad4db65SToby Isaac     PetscInt  i;
199*fad4db65SToby Isaac 
200*fad4db65SToby Isaac     k = PetscMin(k, n - k);
201*fad4db65SToby Isaac     for (i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
202*fad4db65SToby Isaac     *binomial = binom;
203*fad4db65SToby Isaac   }
204*fad4db65SToby Isaac   PetscFunctionReturn(0);
205*fad4db65SToby Isaac }
206*fad4db65SToby Isaac 
207*fad4db65SToby Isaac /*MC
208*fad4db65SToby Isaac    PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps.
209*fad4db65SToby Isaac 
210*fad4db65SToby Isaac    A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation,
211*fad4db65SToby Isaac    by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in
212*fad4db65SToby Isaac    some position j >= i.  We encode this swap as the difference (j - i).  The difference d_i at step i is less than
213*fad4db65SToby Isaac    (n - i).  We encode this sequence of n-1 differences [d_0, ..., d_{n-2}] as the number
214*fad4db65SToby Isaac    (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}.
215*fad4db65SToby Isaac 
216*fad4db65SToby Isaac    Input Arguments:
217*fad4db65SToby Isaac 
218*fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
219*fad4db65SToby Isaac .  k - an integer in [0, n!)
220*fad4db65SToby Isaac -  work - a workspace of n integers
221*fad4db65SToby Isaac 
222*fad4db65SToby Isaac    Output Arguments:
223*fad4db65SToby Isaac 
224*fad4db65SToby Isaac +  perm - the permuted list of the integers [0, ..., n-1]
225*fad4db65SToby Isaac .  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.
226*fad4db65SToby Isaac 
227*fad4db65SToby 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.
228*fad4db65SToby Isaac 
229*fad4db65SToby Isaac    Level: beginner
230*fad4db65SToby Isaac M*/
231*fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd)
2321a989b97SToby Isaac {
2331a989b97SToby Isaac   PetscInt  odd = 0;
2341a989b97SToby Isaac   PetscInt  i;
235*fad4db65SToby Isaac   PetscInt  work[PETSC_FACTORIAL_MAX];
236*fad4db65SToby Isaac   PetscInt *w;
2371a989b97SToby Isaac 
2381a989b97SToby Isaac   PetscFunctionBeginHot;
239*fad4db65SToby 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);
240*fad4db65SToby Isaac   w = &work[n - 2];
2411a989b97SToby Isaac   for (i = 2; i <= n; i++) {
2421a989b97SToby Isaac     *(w--) = k % i;
2431a989b97SToby Isaac     k /= i;
2441a989b97SToby Isaac   }
2451a989b97SToby Isaac   for (i = 0; i < n; i++) perm[i] = i;
2461a989b97SToby Isaac   for (i = 0; i < n - 1; i++) {
2471a989b97SToby Isaac     PetscInt s = work[i];
2481a989b97SToby Isaac     PetscInt swap = perm[i];
2491a989b97SToby Isaac 
2501a989b97SToby Isaac     perm[i] = perm[i + s];
2511a989b97SToby Isaac     perm[i + s] = swap;
2521a989b97SToby Isaac     odd ^= (!!s);
2531a989b97SToby Isaac   }
2541a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
2551a989b97SToby Isaac   PetscFunctionReturn(0);
2561a989b97SToby Isaac }
2571a989b97SToby Isaac 
258*fad4db65SToby Isaac /*MC
259*fad4db65SToby Isaac    PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k).
260*fad4db65SToby Isaac    The encoding is in lexicographic order.
261*fad4db65SToby Isaac 
262*fad4db65SToby Isaac    Input Arguments:
263*fad4db65SToby Isaac 
264*fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
265*fad4db65SToby Isaac .  k - an integer in [0, n]
266*fad4db65SToby Isaac -  j - an index in [0, n choose k)
267*fad4db65SToby Isaac 
268*fad4db65SToby Isaac    Output Arguments:
269*fad4db65SToby Isaac 
270*fad4db65SToby Isaac .  subset - the jth subset of size k of the integers [0, ..., n - 1]
271*fad4db65SToby Isaac 
272*fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
273*fad4db65SToby Isaac 
274*fad4db65SToby Isaac    Level: beginner
275*fad4db65SToby Isaac 
276*fad4db65SToby Isaac .seealso: PetscDTSubsetIndex()
277*fad4db65SToby Isaac M*/
2781a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
2791a989b97SToby Isaac {
2801a989b97SToby Isaac   PetscInt       Nk, i, l;
2811a989b97SToby Isaac   PetscErrorCode ierr;
2821a989b97SToby Isaac 
2831a989b97SToby Isaac   PetscFunctionBeginHot;
284*fad4db65SToby Isaac   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
2851a989b97SToby Isaac   for (i = 0, l = 0; i < n && l < k; i++) {
2861a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
2871a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
2881a989b97SToby Isaac 
2891a989b97SToby Isaac     if (j < Nminuskminus) {
2901a989b97SToby Isaac       subset[l++] = i;
2911a989b97SToby Isaac       Nk = Nminuskminus;
2921a989b97SToby Isaac     } else {
2931a989b97SToby Isaac       j -= Nminuskminus;
2941a989b97SToby Isaac       Nk = Nminusk;
2951a989b97SToby Isaac     }
2961a989b97SToby Isaac   }
2971a989b97SToby Isaac   PetscFunctionReturn(0);
2981a989b97SToby Isaac }
2991a989b97SToby Isaac 
300*fad4db65SToby Isaac /*MC
301*fad4db65SToby 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.
302*fad4db65SToby Isaac 
303*fad4db65SToby Isaac    Input Arguments:
304*fad4db65SToby Isaac 
305*fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
306*fad4db65SToby Isaac .  k - an integer in [0, n]
307*fad4db65SToby Isaac -  subset - an ordered subset of the integers [0, ..., n - 1]
308*fad4db65SToby Isaac 
309*fad4db65SToby Isaac    Output Arguments:
310*fad4db65SToby Isaac 
311*fad4db65SToby Isaac .  index - the rank of the subset in lexicographic order
312*fad4db65SToby Isaac 
313*fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
314*fad4db65SToby Isaac 
315*fad4db65SToby Isaac    Level: beginner
316*fad4db65SToby Isaac 
317*fad4db65SToby Isaac .seealso: PetscDTEnumSubset()
318*fad4db65SToby Isaac M*/
3191a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
3201a989b97SToby Isaac {
3211a989b97SToby Isaac   PetscInt       i, j = 0, l, Nk;
3221a989b97SToby Isaac   PetscErrorCode ierr;
3231a989b97SToby Isaac 
3241a989b97SToby Isaac   PetscFunctionBeginHot;
325*fad4db65SToby Isaac   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
3261a989b97SToby Isaac   for (i = 0, l = 0; i < n && l < k; i++) {
3271a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
3281a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
3291a989b97SToby Isaac 
3301a989b97SToby Isaac     if (subset[l] == i) {
3311a989b97SToby Isaac       l++;
3321a989b97SToby Isaac       Nk = Nminuskminus;
3331a989b97SToby Isaac     } else {
3341a989b97SToby Isaac       j += Nminuskminus;
3351a989b97SToby Isaac       Nk = Nminusk;
3361a989b97SToby Isaac     }
3371a989b97SToby Isaac   }
3381a989b97SToby Isaac   *index = j;
3391a989b97SToby Isaac   PetscFunctionReturn(0);
3401a989b97SToby Isaac }
3411a989b97SToby Isaac 
3421a989b97SToby Isaac 
343*fad4db65SToby Isaac /*MC
344*fad4db65SToby Isaac    PetscDTEnumSubset - Split the integers [0, ..., n - 1] into two complementary ordered subsets, the first of size k and beingthe jth in lexicographic order.
345*fad4db65SToby Isaac 
346*fad4db65SToby Isaac    Input Arguments:
347*fad4db65SToby Isaac 
348*fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
349*fad4db65SToby Isaac .  k - an integer in [0, n]
350*fad4db65SToby Isaac -  j - an index in [0, n choose k)
351*fad4db65SToby Isaac 
352*fad4db65SToby Isaac    Output Arguments:
353*fad4db65SToby Isaac 
354*fad4db65SToby Isaac +  perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set.
355*fad4db65SToby Isaac -  isOdd - if not NULL, return whether the permutation is even or odd.
356*fad4db65SToby Isaac 
357*fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
358*fad4db65SToby Isaac 
359*fad4db65SToby Isaac    Level: beginner
360*fad4db65SToby Isaac 
361*fad4db65SToby Isaac .seealso: PetscDTEnumSubset(), PetscDTSubsetIndex()
362*fad4db65SToby Isaac M*/
363*fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd)
3641a989b97SToby Isaac {
3651a989b97SToby Isaac   PetscInt       i, l, m, *subcomp, Nk;
3661a989b97SToby Isaac   PetscInt       odd;
3671a989b97SToby Isaac   PetscErrorCode ierr;
3681a989b97SToby Isaac 
3691a989b97SToby Isaac   PetscFunctionBeginHot;
370*fad4db65SToby Isaac   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
3711a989b97SToby Isaac   odd = 0;
372*fad4db65SToby Isaac   subcomp = &perm[k];
3731a989b97SToby Isaac   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
3741a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
3751a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
3761a989b97SToby Isaac 
3771a989b97SToby Isaac     if (j < Nminuskminus) {
378*fad4db65SToby Isaac       perm[l++] = i;
3791a989b97SToby Isaac       Nk = Nminuskminus;
3801a989b97SToby Isaac     } else {
3811a989b97SToby Isaac       subcomp[m++] = i;
3821a989b97SToby Isaac       j -= Nminuskminus;
3831a989b97SToby Isaac       odd ^= ((k - l) & 1);
3841a989b97SToby Isaac       Nk = Nminusk;
3851a989b97SToby Isaac     }
3861a989b97SToby Isaac   }
3871a989b97SToby Isaac   for (; i < n; i++) {
3881a989b97SToby Isaac     subcomp[m++] = i;
3891a989b97SToby Isaac   }
3901a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
3911a989b97SToby Isaac   PetscFunctionReturn(0);
3921a989b97SToby Isaac }
3931a989b97SToby Isaac 
39437045ce4SJed Brown #endif
395