xref: /petsc/include/petscdt.h (revision e2ab39cceefafbe1e2ff2e70de80588e875d332a)
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 
92cd22861SMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCQUADRATURE_CLASSID;
102cd22861SMatthew G. Knepley 
1121454ff5SMatthew G. Knepley /*S
1221454ff5SMatthew G. Knepley   PetscQuadrature - Quadrature rule for integration.
1321454ff5SMatthew G. Knepley 
14329bbf4eSMatthew G. Knepley   Level: beginner
1521454ff5SMatthew G. Knepley 
1621454ff5SMatthew G. Knepley .seealso:  PetscQuadratureCreate(), PetscQuadratureDestroy()
1721454ff5SMatthew G. Knepley S*/
1821454ff5SMatthew G. Knepley typedef struct _p_PetscQuadrature *PetscQuadrature;
1921454ff5SMatthew G. Knepley 
208272889dSSatish Balay /*E
21916e780bShannah_mairs   PetscGaussLobattoLegendreCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights
228272889dSSatish Balay 
238272889dSSatish Balay   Level: intermediate
248272889dSSatish Balay 
25f2e8fe4dShannah_mairs $  PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA - compute the nodes via linear algebra
26d410ae54Shannah_mairs $  PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON - compute the nodes by solving a nonlinear equation with Newton's method
278272889dSSatish Balay 
288272889dSSatish Balay E*/
29f2e8fe4dShannah_mairs typedef enum {PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON} PetscGaussLobattoLegendreCreateType;
308272889dSSatish Balay 
3121454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *);
32c9638911SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *);
33bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt*);
34bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt);
35a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt*);
36a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt);
37a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt*, PetscInt*, PetscInt*, const PetscReal *[], const PetscReal *[]);
38a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal [], const PetscReal []);
3921454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer);
4021454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *);
41a0845e3aSMatthew G. Knepley 
4289710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *);
4389710940SMatthew G. Knepley 
44907761f8SToby Isaac PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *);
45907761f8SToby Isaac 
4637045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*);
4737045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt,PetscReal,PetscReal,PetscReal*,PetscReal*);
48916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt,PetscGaussLobattoLegendreCreateType,PetscReal*,PetscReal*);
49194825f6SJed Brown PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt,PetscInt,const PetscReal*,PetscInt,const PetscReal*,PetscReal*);
50a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
51a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
5237045ce4SJed Brown 
53b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
54b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *);
55d525116cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *);
56b3c0f97bSTom Klotz 
57916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *);
58916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
59916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
60916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
61916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
62916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
63916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
64916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
65916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
66916e780bShannah_mairs 
671a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
681a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
691a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
701a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
711a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
721a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
731a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
74dda711d0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
751a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
761a989b97SToby Isaac 
77fad4db65SToby Isaac #if defined(PETSC_USE_64BIT_INDICES)
78fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 20
79fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX  61
80fad4db65SToby Isaac #else
81fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 12
82fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX  29
83fad4db65SToby Isaac #endif
84fad4db65SToby Isaac 
85fad4db65SToby Isaac /*MC
86fad4db65SToby Isaac    PetscDTFactorial - Approximate n! as a real number
87fad4db65SToby Isaac 
88fad4db65SToby Isaac    Input Arguments:
89fad4db65SToby Isaac .  n - a non-negative integer
90fad4db65SToby Isaac 
9128222859SToby Isaac    Output Arguments:
92fad4db65SToby Isaac .  factorial - n!
93fad4db65SToby Isaac 
94fad4db65SToby Isaac    Level: beginner
95fad4db65SToby Isaac M*/
96fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial)
97fad4db65SToby Isaac {
98fad4db65SToby Isaac   PetscReal f = 1.0;
99fad4db65SToby Isaac   PetscInt  i;
100fad4db65SToby Isaac 
101fad4db65SToby Isaac   PetscFunctionBegin;
102*e2ab39ccSLisandro Dalcin   *factorial = -1.0;
10328222859SToby Isaac   if (n < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %D\n", n);
104*e2ab39ccSLisandro Dalcin   for (i = 1; i < n+1; ++i) f *= (PetscReal)i;
105fad4db65SToby Isaac   *factorial = f;
106fad4db65SToby Isaac   PetscFunctionReturn(0);
107fad4db65SToby Isaac }
108fad4db65SToby Isaac 
109fad4db65SToby Isaac /*MC
110fad4db65SToby Isaac    PetscDTFactorialInt - Compute n! as an integer
111fad4db65SToby Isaac 
112fad4db65SToby Isaac    Input Arguments:
113fad4db65SToby Isaac .  n - a non-negative integer
114fad4db65SToby Isaac 
11528222859SToby Isaac    Output Arguments:
116fad4db65SToby Isaac .  factorial - n!
117fad4db65SToby Isaac 
118fad4db65SToby Isaac    Level: beginner
119fad4db65SToby Isaac 
120fad4db65SToby 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.
121fad4db65SToby Isaac M*/
122fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial)
123fad4db65SToby Isaac {
124fad4db65SToby Isaac   PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600};
125fad4db65SToby Isaac 
12628222859SToby Isaac   PetscFunctionBegin;
12728222859SToby Isaac   *factorial = -1;
128fad4db65SToby 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);
129fad4db65SToby Isaac   if (n <= 12) {
130fad4db65SToby Isaac     *factorial = facLookup[n];
131fad4db65SToby Isaac   } else {
132fad4db65SToby Isaac     PetscInt f = facLookup[12];
133fad4db65SToby Isaac     PetscInt i;
134fad4db65SToby Isaac 
135fad4db65SToby Isaac     for (i = 13; i < n+1; ++i) f *= i;
136fad4db65SToby Isaac     *factorial = f;
137fad4db65SToby Isaac   }
138fad4db65SToby Isaac   PetscFunctionReturn(0);
139fad4db65SToby Isaac }
140fad4db65SToby Isaac 
141fad4db65SToby Isaac /*MC
142fad4db65SToby Isaac    PetscDTBinomial - Approximate the binomial coefficient "n choose k"
143fad4db65SToby Isaac 
144fad4db65SToby Isaac    Input Arguments:
145fad4db65SToby Isaac +  n - a non-negative integer
146fad4db65SToby Isaac -  k - an integer between 0 and n, inclusive
147fad4db65SToby Isaac 
14828222859SToby Isaac    Output Arguments:
149fad4db65SToby Isaac .  binomial - approximation of the binomial coefficient n choose k
150fad4db65SToby Isaac 
151fad4db65SToby Isaac    Level: beginner
152fad4db65SToby Isaac M*/
153fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
1541a989b97SToby Isaac {
1551a989b97SToby Isaac   PetscFunctionBeginHot;
156*e2ab39ccSLisandro Dalcin   *binomial = -1.0;
157fad4db65SToby 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);
1581a989b97SToby Isaac   if (n <= 3) {
1591a989b97SToby Isaac     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};
1601a989b97SToby Isaac 
161*e2ab39ccSLisandro Dalcin     *binomial = (PetscReal)binomLookup[n][k];
1621a989b97SToby Isaac   } else {
163*e2ab39ccSLisandro Dalcin     PetscReal binom = 1.0;
1641a989b97SToby Isaac     PetscInt  i;
1651a989b97SToby Isaac 
1661a989b97SToby Isaac     k = PetscMin(k, n - k);
167*e2ab39ccSLisandro Dalcin     for (i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1);
1681a989b97SToby Isaac     *binomial = binom;
1691a989b97SToby Isaac   }
1701a989b97SToby Isaac   PetscFunctionReturn(0);
1711a989b97SToby Isaac }
1721a989b97SToby Isaac 
173fad4db65SToby Isaac /*MC
174fad4db65SToby Isaac    PetscDTBinomialInt - Compute the binomial coefficient "n choose k"
175fad4db65SToby Isaac 
176fad4db65SToby Isaac    Input Arguments:
177fad4db65SToby Isaac +  n - a non-negative integer
178fad4db65SToby Isaac -  k - an integer between 0 and n, inclusive
179fad4db65SToby Isaac 
18028222859SToby Isaac    Output Arguments:
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 {
18928222859SToby Isaac   PetscInt bin;
19028222859SToby Isaac 
19128222859SToby Isaac   PetscFunctionBegin;
19228222859SToby Isaac   *binomial = -1;
193fad4db65SToby 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);
194fad4db65SToby 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);
195fad4db65SToby Isaac   if (n <= 3) {
196fad4db65SToby Isaac     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};
197fad4db65SToby Isaac 
19828222859SToby Isaac     bin = binomLookup[n][k];
199fad4db65SToby Isaac   } else {
200fad4db65SToby Isaac     PetscInt  binom = 1;
201fad4db65SToby Isaac     PetscInt  i;
202fad4db65SToby Isaac 
203fad4db65SToby Isaac     k = PetscMin(k, n - k);
204fad4db65SToby Isaac     for (i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
20528222859SToby Isaac     bin = binom;
206fad4db65SToby Isaac   }
20728222859SToby Isaac   *binomial = bin;
208fad4db65SToby Isaac   PetscFunctionReturn(0);
209fad4db65SToby Isaac }
210fad4db65SToby Isaac 
211fad4db65SToby Isaac /*MC
212fad4db65SToby Isaac    PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps.
213fad4db65SToby Isaac 
214fad4db65SToby Isaac    A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation,
215fad4db65SToby Isaac    by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in
21628222859SToby Isaac    some position j >= i.  This swap is encoded as the difference (j - i).  The difference d_i at step i is less than
21728222859SToby Isaac    (n - i).  This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number
218fad4db65SToby Isaac    (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}.
219fad4db65SToby Isaac 
220fad4db65SToby Isaac    Input Arguments:
221fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
2228cd1e013SToby Isaac -  k - an integer in [0, n!)
223fad4db65SToby Isaac 
224fad4db65SToby Isaac    Output Arguments:
225fad4db65SToby Isaac +  perm - the permuted list of the integers [0, ..., n-1]
2268cd1e013SToby Isaac -  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.
227fad4db65SToby Isaac 
228fad4db65SToby 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.
229fad4db65SToby Isaac 
230fad4db65SToby Isaac    Level: beginner
231fad4db65SToby Isaac M*/
232fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd)
2331a989b97SToby Isaac {
2341a989b97SToby Isaac   PetscInt  odd = 0;
2351a989b97SToby Isaac   PetscInt  i;
236fad4db65SToby Isaac   PetscInt  work[PETSC_FACTORIAL_MAX];
237fad4db65SToby Isaac   PetscInt *w;
2381a989b97SToby Isaac 
23928222859SToby Isaac   PetscFunctionBegin;
24028222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
241fad4db65SToby 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);
242fad4db65SToby Isaac   w = &work[n - 2];
2431a989b97SToby Isaac   for (i = 2; i <= n; i++) {
2441a989b97SToby Isaac     *(w--) = k % i;
2451a989b97SToby Isaac     k /= i;
2461a989b97SToby Isaac   }
2471a989b97SToby Isaac   for (i = 0; i < n; i++) perm[i] = i;
2481a989b97SToby Isaac   for (i = 0; i < n - 1; i++) {
2491a989b97SToby Isaac     PetscInt s = work[i];
2501a989b97SToby Isaac     PetscInt swap = perm[i];
2511a989b97SToby Isaac 
2521a989b97SToby Isaac     perm[i] = perm[i + s];
2531a989b97SToby Isaac     perm[i + s] = swap;
2541a989b97SToby Isaac     odd ^= (!!s);
2551a989b97SToby Isaac   }
2561a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
2571a989b97SToby Isaac   PetscFunctionReturn(0);
2581a989b97SToby Isaac }
2591a989b97SToby Isaac 
260fad4db65SToby Isaac /*MC
2618cd1e013SToby Isaac    PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!).  This inverts PetscDTEnumPerm.
2628cd1e013SToby Isaac 
2638cd1e013SToby Isaac    Input Arguments:
2648cd1e013SToby Isaac +  n - a non-negative integer (see note about limits below)
2658cd1e013SToby Isaac -  perm - the permuted list of the integers [0, ..., n-1]
2668cd1e013SToby Isaac 
2678cd1e013SToby Isaac    Output Arguments:
2688cd1e013SToby Isaac +  k - an integer in [0, n!)
2698cd1e013SToby Isaac .  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.
2708cd1e013SToby Isaac 
2718cd1e013SToby 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.
2728cd1e013SToby Isaac 
2738cd1e013SToby Isaac    Level: beginner
2748cd1e013SToby Isaac M*/
2758cd1e013SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd)
2768cd1e013SToby Isaac {
2778cd1e013SToby Isaac   PetscInt  odd = 0;
2788cd1e013SToby Isaac   PetscInt  i, idx;
2798cd1e013SToby Isaac   PetscInt  work[PETSC_FACTORIAL_MAX];
2808cd1e013SToby Isaac   PetscInt  iwork[PETSC_FACTORIAL_MAX];
2818cd1e013SToby Isaac 
2828cd1e013SToby Isaac   PetscFunctionBeginHot;
28328222859SToby Isaac   *k = -1;
28428222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
2858cd1e013SToby 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);
2868cd1e013SToby Isaac   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
2878cd1e013SToby Isaac   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
2888cd1e013SToby Isaac   for (idx = 0, i = 0; i < n - 1; i++) {
2898cd1e013SToby Isaac     PetscInt j = perm[i];
2908cd1e013SToby Isaac     PetscInt icur = work[i];
2918cd1e013SToby Isaac     PetscInt jloc = iwork[j];
2928cd1e013SToby Isaac     PetscInt diff = jloc - i;
2938cd1e013SToby Isaac 
2948cd1e013SToby Isaac     idx = idx * (n - i) + diff;
2958cd1e013SToby Isaac     /* swap (i, jloc) */
2968cd1e013SToby Isaac     work[i] = j;
2978cd1e013SToby Isaac     work[jloc] = icur;
2988cd1e013SToby Isaac     iwork[j] = i;
2998cd1e013SToby Isaac     iwork[icur] = jloc;
3008cd1e013SToby Isaac     odd ^= (!!diff);
3018cd1e013SToby Isaac   }
3028cd1e013SToby Isaac   *k = idx;
3038cd1e013SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
3048cd1e013SToby Isaac   PetscFunctionReturn(0);
3058cd1e013SToby Isaac }
3068cd1e013SToby Isaac 
3078cd1e013SToby Isaac /*MC
308fad4db65SToby Isaac    PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k).
309fad4db65SToby Isaac    The encoding is in lexicographic order.
310fad4db65SToby Isaac 
311fad4db65SToby Isaac    Input Arguments:
312fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
313fad4db65SToby Isaac .  k - an integer in [0, n]
314fad4db65SToby Isaac -  j - an index in [0, n choose k)
315fad4db65SToby Isaac 
316fad4db65SToby Isaac    Output Arguments:
317fad4db65SToby Isaac .  subset - the jth subset of size k of the integers [0, ..., n - 1]
318fad4db65SToby Isaac 
319fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
320fad4db65SToby Isaac 
321fad4db65SToby Isaac    Level: beginner
322fad4db65SToby Isaac 
323fad4db65SToby Isaac .seealso: PetscDTSubsetIndex()
324fad4db65SToby Isaac M*/
3251a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
3261a989b97SToby Isaac {
3271a989b97SToby Isaac   PetscInt       Nk, i, l;
3281a989b97SToby Isaac   PetscErrorCode ierr;
3291a989b97SToby Isaac 
3301a989b97SToby Isaac   PetscFunctionBeginHot;
331fad4db65SToby Isaac   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
3321a989b97SToby Isaac   for (i = 0, l = 0; i < n && l < k; i++) {
3331a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
3341a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
3351a989b97SToby Isaac 
3361a989b97SToby Isaac     if (j < Nminuskminus) {
3371a989b97SToby Isaac       subset[l++] = i;
3381a989b97SToby Isaac       Nk = Nminuskminus;
3391a989b97SToby Isaac     } else {
3401a989b97SToby Isaac       j -= Nminuskminus;
3411a989b97SToby Isaac       Nk = Nminusk;
3421a989b97SToby Isaac     }
3431a989b97SToby Isaac   }
3441a989b97SToby Isaac   PetscFunctionReturn(0);
3451a989b97SToby Isaac }
3461a989b97SToby Isaac 
347fad4db65SToby Isaac /*MC
348fad4db65SToby 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.
349fad4db65SToby Isaac 
350fad4db65SToby Isaac    Input Arguments:
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 .  index - the rank of the subset in lexicographic order
357fad4db65SToby Isaac 
358fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
359fad4db65SToby Isaac 
360fad4db65SToby Isaac    Level: beginner
361fad4db65SToby Isaac 
362fad4db65SToby Isaac .seealso: PetscDTEnumSubset()
363fad4db65SToby Isaac M*/
3641a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
3651a989b97SToby Isaac {
3661a989b97SToby Isaac   PetscInt       i, j = 0, l, Nk;
3671a989b97SToby Isaac   PetscErrorCode ierr;
3681a989b97SToby Isaac 
36928222859SToby Isaac   PetscFunctionBegin;
37028222859SToby Isaac   *index = -1;
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 
388fad4db65SToby Isaac /*MC
38928222859SToby 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.
390fad4db65SToby Isaac 
391fad4db65SToby Isaac    Input Arguments:
392fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
393fad4db65SToby Isaac .  k - an integer in [0, n]
394fad4db65SToby Isaac -  j - an index in [0, n choose k)
395fad4db65SToby Isaac 
396fad4db65SToby Isaac    Output Arguments:
397fad4db65SToby Isaac +  perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set.
39828222859SToby Isaac -  isOdd - if not NULL, return whether perm is an even or odd permutation.
399fad4db65SToby Isaac 
400fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
401fad4db65SToby Isaac 
402fad4db65SToby Isaac    Level: beginner
403fad4db65SToby Isaac 
404fad4db65SToby Isaac .seealso: PetscDTEnumSubset(), PetscDTSubsetIndex()
405fad4db65SToby Isaac M*/
406fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd)
4071a989b97SToby Isaac {
4081a989b97SToby Isaac   PetscInt       i, l, m, *subcomp, Nk;
4091a989b97SToby Isaac   PetscInt       odd;
4101a989b97SToby Isaac   PetscErrorCode ierr;
4111a989b97SToby Isaac 
41228222859SToby Isaac   PetscFunctionBegin;
41328222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
414fad4db65SToby Isaac   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
4151a989b97SToby Isaac   odd = 0;
416fad4db65SToby Isaac   subcomp = &perm[k];
4171a989b97SToby Isaac   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
4181a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
4191a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
4201a989b97SToby Isaac 
4211a989b97SToby Isaac     if (j < Nminuskminus) {
422fad4db65SToby Isaac       perm[l++] = i;
4231a989b97SToby Isaac       Nk = Nminuskminus;
4241a989b97SToby Isaac     } else {
4251a989b97SToby Isaac       subcomp[m++] = i;
4261a989b97SToby Isaac       j -= Nminuskminus;
4271a989b97SToby Isaac       odd ^= ((k - l) & 1);
4281a989b97SToby Isaac       Nk = Nminusk;
4291a989b97SToby Isaac     }
4301a989b97SToby Isaac   }
4311a989b97SToby Isaac   for (; i < n; i++) {
4321a989b97SToby Isaac     subcomp[m++] = i;
4331a989b97SToby Isaac   }
4341a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
4351a989b97SToby Isaac   PetscFunctionReturn(0);
4361a989b97SToby Isaac }
4371a989b97SToby Isaac 
438ef0bb6c7SMatthew G. Knepley struct _p_PetscTabulation {
439ef0bb6c7SMatthew G. Knepley   PetscInt    K;    /* Indicates a k-jet, namely tabulated derviatives up to order k */
440ef0bb6c7SMatthew G. Knepley   PetscInt    Nr;   /* THe number of tabulation replicas (often 1) */
441ef0bb6c7SMatthew G. Knepley   PetscInt    Np;   /* The number of tabulation points in a replica */
442ef0bb6c7SMatthew G. Knepley   PetscInt    Nb;   /* The number of functions tabulated */
443ef0bb6c7SMatthew G. Knepley   PetscInt    Nc;   /* The number of function components */
444ef0bb6c7SMatthew G. Knepley   PetscInt    cdim; /* The coordinate dimension */
445ef0bb6c7SMatthew G. Knepley   PetscReal **T;    /* The tabulation T[K] of functions and their derivatives
446ef0bb6c7SMatthew G. Knepley                        T[0] = B[Nr*Np][Nb][Nc]:             The basis function values at quadrature points
447ef0bb6c7SMatthew G. Knepley                        T[1] = D[Nr*Np][Nb][Nc][cdim]:       The basis function derivatives at quadrature points
448ef0bb6c7SMatthew G. Knepley                        T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */
449ef0bb6c7SMatthew G. Knepley };
450ef0bb6c7SMatthew G. Knepley typedef struct _p_PetscTabulation *PetscTabulation;
451ef0bb6c7SMatthew G. Knepley 
45237045ce4SJed Brown #endif
453