xref: /petsc/include/petscdt.h (revision d4afb7202e2a84f4aed53e5b2470b60516492222)
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 
31*d4afb720SToby Isaac /*E
32*d4afb720SToby Isaac   PetscDTNodeType - A description of strategies for generating nodes (both
33*d4afb720SToby Isaac   quadrature nodes and nodes for Lagrange polynomials)
34*d4afb720SToby Isaac 
35*d4afb720SToby Isaac   Level: intermediate
36*d4afb720SToby Isaac 
37*d4afb720SToby Isaac $  PETSCDTNODES_DEFAULT - Nodes chosen by PETSc
38*d4afb720SToby Isaac $  PETSCDTNODES_GAUSSJACOBI - Nodes at either Gauss-Jacobi or Gauss-Lobatto-Jacobi quadrature points
39*d4afb720SToby Isaac $  PETSCDTNODES_EQUISPACED - Nodes equispaced either including the endpoints or excluding them
40*d4afb720SToby Isaac $  PETSCDTNODES_TANHSINH - Nodes at Tanh-Sinh quadrature points
41*d4afb720SToby Isaac 
42*d4afb720SToby Isaac   Note: a PetscDTNodeType can be paired with a PetscBool to indicate whether
43*d4afb720SToby Isaac   the nodes include endpoints or not, and in the case of PETSCDT_GAUSSJACOBI
44*d4afb720SToby Isaac   with exponents for the weight function.
45*d4afb720SToby Isaac 
46*d4afb720SToby Isaac E*/
47*d4afb720SToby Isaac typedef enum {PETSCDTNODES_DEFAULT=-1, PETSCDTNODES_GAUSSJACOBI, PETSCDTNODES_EQUISPACED, PETSCDTNODES_TANHSINH} PetscDTNodeType;
48*d4afb720SToby Isaac 
49*d4afb720SToby Isaac PETSC_EXTERN const char *const PetscDTNodeTypes[];
50*d4afb720SToby Isaac 
5121454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *);
52c9638911SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *);
53bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt*);
54bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt);
55a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt*);
56a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt);
57a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt*, PetscInt*, PetscInt*, const PetscReal *[], const PetscReal *[]);
58a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal [], const PetscReal []);
5921454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer);
6021454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *);
61a0845e3aSMatthew G. Knepley 
6289710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *);
6389710940SMatthew G. Knepley 
64907761f8SToby Isaac PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *);
65907761f8SToby Isaac 
6637045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*);
6794e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt,PetscReal,PetscReal,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*);
6837045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt,PetscReal,PetscReal,PetscReal*,PetscReal*);
6994e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*,PetscReal*);
7094e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*,PetscReal*);
71916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt,PetscGaussLobattoLegendreCreateType,PetscReal*,PetscReal*);
72194825f6SJed Brown PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt,PetscInt,const PetscReal*,PetscInt,const PetscReal*,PetscReal*);
73a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
74e6a796c3SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
7537045ce4SJed Brown 
76b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
77b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *);
78d525116cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *);
79b3c0f97bSTom Klotz 
80916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *);
81916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
82916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
83916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
84916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
85916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
86916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
87916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
88916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
89916e780bShannah_mairs 
901a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
911a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
921a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
931a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
941a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
951a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
961a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
97dda711d0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
981a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
991a989b97SToby Isaac 
100*d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt,PetscInt,const PetscInt[],PetscInt*);
101*d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt,PetscInt,PetscInt,PetscInt[]);
102*d4afb720SToby Isaac 
103fad4db65SToby Isaac #if defined(PETSC_USE_64BIT_INDICES)
104fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 20
105fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX  61
106fad4db65SToby Isaac #else
107fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 12
108fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX  29
109fad4db65SToby Isaac #endif
110fad4db65SToby Isaac 
111fad4db65SToby Isaac /*MC
112fad4db65SToby Isaac    PetscDTFactorial - Approximate n! as a real number
113fad4db65SToby Isaac 
114fad4db65SToby Isaac    Input Arguments:
115fad4db65SToby Isaac .  n - a non-negative integer
116fad4db65SToby Isaac 
11728222859SToby Isaac    Output Arguments:
118fad4db65SToby Isaac .  factorial - n!
119fad4db65SToby Isaac 
120fad4db65SToby Isaac    Level: beginner
121fad4db65SToby Isaac M*/
122fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial)
123fad4db65SToby Isaac {
124fad4db65SToby Isaac   PetscReal f = 1.0;
125fad4db65SToby Isaac   PetscInt  i;
126fad4db65SToby Isaac 
127fad4db65SToby Isaac   PetscFunctionBegin;
128e2ab39ccSLisandro Dalcin   *factorial = -1.0;
12928222859SToby Isaac   if (n < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %D\n", n);
130e2ab39ccSLisandro Dalcin   for (i = 1; i < n+1; ++i) f *= (PetscReal)i;
131fad4db65SToby Isaac   *factorial = f;
132fad4db65SToby Isaac   PetscFunctionReturn(0);
133fad4db65SToby Isaac }
134fad4db65SToby Isaac 
135fad4db65SToby Isaac /*MC
136fad4db65SToby Isaac    PetscDTFactorialInt - Compute n! as an integer
137fad4db65SToby Isaac 
138fad4db65SToby Isaac    Input Arguments:
139fad4db65SToby Isaac .  n - a non-negative integer
140fad4db65SToby Isaac 
14128222859SToby Isaac    Output Arguments:
142fad4db65SToby Isaac .  factorial - n!
143fad4db65SToby Isaac 
144fad4db65SToby Isaac    Level: beginner
145fad4db65SToby Isaac 
146fad4db65SToby 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.
147fad4db65SToby Isaac M*/
148fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial)
149fad4db65SToby Isaac {
150fad4db65SToby Isaac   PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600};
151fad4db65SToby Isaac 
15228222859SToby Isaac   PetscFunctionBegin;
15328222859SToby Isaac   *factorial = -1;
154fad4db65SToby 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);
155fad4db65SToby Isaac   if (n <= 12) {
156fad4db65SToby Isaac     *factorial = facLookup[n];
157fad4db65SToby Isaac   } else {
158fad4db65SToby Isaac     PetscInt f = facLookup[12];
159fad4db65SToby Isaac     PetscInt i;
160fad4db65SToby Isaac 
161fad4db65SToby Isaac     for (i = 13; i < n+1; ++i) f *= i;
162fad4db65SToby Isaac     *factorial = f;
163fad4db65SToby Isaac   }
164fad4db65SToby Isaac   PetscFunctionReturn(0);
165fad4db65SToby Isaac }
166fad4db65SToby Isaac 
167fad4db65SToby Isaac /*MC
168fad4db65SToby Isaac    PetscDTBinomial - Approximate the binomial coefficient "n choose k"
169fad4db65SToby Isaac 
170fad4db65SToby Isaac    Input Arguments:
171fad4db65SToby Isaac +  n - a non-negative integer
172fad4db65SToby Isaac -  k - an integer between 0 and n, inclusive
173fad4db65SToby Isaac 
17428222859SToby Isaac    Output Arguments:
175fad4db65SToby Isaac .  binomial - approximation of the binomial coefficient n choose k
176fad4db65SToby Isaac 
177fad4db65SToby Isaac    Level: beginner
178fad4db65SToby Isaac M*/
179fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
1801a989b97SToby Isaac {
1811a989b97SToby Isaac   PetscFunctionBeginHot;
182e2ab39ccSLisandro Dalcin   *binomial = -1.0;
183fad4db65SToby 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);
1841a989b97SToby Isaac   if (n <= 3) {
1851a989b97SToby Isaac     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};
1861a989b97SToby Isaac 
187e2ab39ccSLisandro Dalcin     *binomial = (PetscReal)binomLookup[n][k];
1881a989b97SToby Isaac   } else {
189e2ab39ccSLisandro Dalcin     PetscReal binom = 1.0;
1901a989b97SToby Isaac     PetscInt  i;
1911a989b97SToby Isaac 
1921a989b97SToby Isaac     k = PetscMin(k, n - k);
193e2ab39ccSLisandro Dalcin     for (i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1);
1941a989b97SToby Isaac     *binomial = binom;
1951a989b97SToby Isaac   }
1961a989b97SToby Isaac   PetscFunctionReturn(0);
1971a989b97SToby Isaac }
1981a989b97SToby Isaac 
199fad4db65SToby Isaac /*MC
200fad4db65SToby Isaac    PetscDTBinomialInt - Compute the binomial coefficient "n choose k"
201fad4db65SToby Isaac 
202fad4db65SToby Isaac    Input Arguments:
203fad4db65SToby Isaac +  n - a non-negative integer
204fad4db65SToby Isaac -  k - an integer between 0 and n, inclusive
205fad4db65SToby Isaac 
20628222859SToby Isaac    Output Arguments:
207fad4db65SToby Isaac .  binomial - the binomial coefficient n choose k
208fad4db65SToby Isaac 
209fad4db65SToby Isaac    Note: this is limited by integers that can be represented by PetscInt
210fad4db65SToby Isaac 
211fad4db65SToby Isaac    Level: beginner
212fad4db65SToby Isaac M*/
213fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial)
214fad4db65SToby Isaac {
21528222859SToby Isaac   PetscInt bin;
21628222859SToby Isaac 
21728222859SToby Isaac   PetscFunctionBegin;
21828222859SToby Isaac   *binomial = -1;
219fad4db65SToby 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);
220fad4db65SToby 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);
221fad4db65SToby Isaac   if (n <= 3) {
222fad4db65SToby Isaac     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};
223fad4db65SToby Isaac 
22428222859SToby Isaac     bin = binomLookup[n][k];
225fad4db65SToby Isaac   } else {
226fad4db65SToby Isaac     PetscInt  binom = 1;
227fad4db65SToby Isaac     PetscInt  i;
228fad4db65SToby Isaac 
229fad4db65SToby Isaac     k = PetscMin(k, n - k);
230fad4db65SToby Isaac     for (i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
23128222859SToby Isaac     bin = binom;
232fad4db65SToby Isaac   }
23328222859SToby Isaac   *binomial = bin;
234fad4db65SToby Isaac   PetscFunctionReturn(0);
235fad4db65SToby Isaac }
236fad4db65SToby Isaac 
237fad4db65SToby Isaac /*MC
238fad4db65SToby Isaac    PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps.
239fad4db65SToby Isaac 
240fad4db65SToby Isaac    A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation,
241fad4db65SToby Isaac    by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in
24228222859SToby Isaac    some position j >= i.  This swap is encoded as the difference (j - i).  The difference d_i at step i is less than
24328222859SToby Isaac    (n - i).  This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number
244fad4db65SToby Isaac    (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}.
245fad4db65SToby Isaac 
246fad4db65SToby Isaac    Input Arguments:
247fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
2488cd1e013SToby Isaac -  k - an integer in [0, n!)
249fad4db65SToby Isaac 
250fad4db65SToby Isaac    Output Arguments:
251fad4db65SToby Isaac +  perm - the permuted list of the integers [0, ..., n-1]
2528cd1e013SToby Isaac -  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.
253fad4db65SToby Isaac 
254fad4db65SToby 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.
255fad4db65SToby Isaac 
256fad4db65SToby Isaac    Level: beginner
257fad4db65SToby Isaac M*/
258fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd)
2591a989b97SToby Isaac {
2601a989b97SToby Isaac   PetscInt  odd = 0;
2611a989b97SToby Isaac   PetscInt  i;
262fad4db65SToby Isaac   PetscInt  work[PETSC_FACTORIAL_MAX];
263fad4db65SToby Isaac   PetscInt *w;
2641a989b97SToby Isaac 
26528222859SToby Isaac   PetscFunctionBegin;
26628222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
267fad4db65SToby 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);
268fad4db65SToby Isaac   w = &work[n - 2];
2691a989b97SToby Isaac   for (i = 2; i <= n; i++) {
2701a989b97SToby Isaac     *(w--) = k % i;
2711a989b97SToby Isaac     k /= i;
2721a989b97SToby Isaac   }
2731a989b97SToby Isaac   for (i = 0; i < n; i++) perm[i] = i;
2741a989b97SToby Isaac   for (i = 0; i < n - 1; i++) {
2751a989b97SToby Isaac     PetscInt s = work[i];
2761a989b97SToby Isaac     PetscInt swap = perm[i];
2771a989b97SToby Isaac 
2781a989b97SToby Isaac     perm[i] = perm[i + s];
2791a989b97SToby Isaac     perm[i + s] = swap;
2801a989b97SToby Isaac     odd ^= (!!s);
2811a989b97SToby Isaac   }
2821a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
2831a989b97SToby Isaac   PetscFunctionReturn(0);
2841a989b97SToby Isaac }
2851a989b97SToby Isaac 
286fad4db65SToby Isaac /*MC
2878cd1e013SToby Isaac    PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!).  This inverts PetscDTEnumPerm.
2888cd1e013SToby Isaac 
2898cd1e013SToby Isaac    Input Arguments:
2908cd1e013SToby Isaac +  n - a non-negative integer (see note about limits below)
2918cd1e013SToby Isaac -  perm - the permuted list of the integers [0, ..., n-1]
2928cd1e013SToby Isaac 
2938cd1e013SToby Isaac    Output Arguments:
2948cd1e013SToby Isaac +  k - an integer in [0, n!)
295f0fc11ceSJed Brown -  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.
2968cd1e013SToby Isaac 
2978cd1e013SToby 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.
2988cd1e013SToby Isaac 
2998cd1e013SToby Isaac    Level: beginner
3008cd1e013SToby Isaac M*/
3018cd1e013SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd)
3028cd1e013SToby Isaac {
3038cd1e013SToby Isaac   PetscInt  odd = 0;
3048cd1e013SToby Isaac   PetscInt  i, idx;
3058cd1e013SToby Isaac   PetscInt  work[PETSC_FACTORIAL_MAX];
3068cd1e013SToby Isaac   PetscInt  iwork[PETSC_FACTORIAL_MAX];
3078cd1e013SToby Isaac 
3088cd1e013SToby Isaac   PetscFunctionBeginHot;
30928222859SToby Isaac   *k = -1;
31028222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
3118cd1e013SToby 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);
3128cd1e013SToby Isaac   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
3138cd1e013SToby Isaac   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
3148cd1e013SToby Isaac   for (idx = 0, i = 0; i < n - 1; i++) {
3158cd1e013SToby Isaac     PetscInt j = perm[i];
3168cd1e013SToby Isaac     PetscInt icur = work[i];
3178cd1e013SToby Isaac     PetscInt jloc = iwork[j];
3188cd1e013SToby Isaac     PetscInt diff = jloc - i;
3198cd1e013SToby Isaac 
3208cd1e013SToby Isaac     idx = idx * (n - i) + diff;
3218cd1e013SToby Isaac     /* swap (i, jloc) */
3228cd1e013SToby Isaac     work[i] = j;
3238cd1e013SToby Isaac     work[jloc] = icur;
3248cd1e013SToby Isaac     iwork[j] = i;
3258cd1e013SToby Isaac     iwork[icur] = jloc;
3268cd1e013SToby Isaac     odd ^= (!!diff);
3278cd1e013SToby Isaac   }
3288cd1e013SToby Isaac   *k = idx;
3298cd1e013SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
3308cd1e013SToby Isaac   PetscFunctionReturn(0);
3318cd1e013SToby Isaac }
3328cd1e013SToby Isaac 
3338cd1e013SToby Isaac /*MC
334fad4db65SToby Isaac    PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k).
335fad4db65SToby Isaac    The encoding is in lexicographic order.
336fad4db65SToby Isaac 
337fad4db65SToby Isaac    Input Arguments:
338fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
339fad4db65SToby Isaac .  k - an integer in [0, n]
340fad4db65SToby Isaac -  j - an index in [0, n choose k)
341fad4db65SToby Isaac 
342fad4db65SToby Isaac    Output Arguments:
343fad4db65SToby Isaac .  subset - the jth subset of size k of the integers [0, ..., n - 1]
344fad4db65SToby Isaac 
345fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
346fad4db65SToby Isaac 
347fad4db65SToby Isaac    Level: beginner
348fad4db65SToby Isaac 
349fad4db65SToby Isaac .seealso: PetscDTSubsetIndex()
350fad4db65SToby Isaac M*/
3511a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
3521a989b97SToby Isaac {
3531a989b97SToby Isaac   PetscInt       Nk, i, l;
3541a989b97SToby Isaac   PetscErrorCode ierr;
3551a989b97SToby Isaac 
3561a989b97SToby Isaac   PetscFunctionBeginHot;
357fad4db65SToby Isaac   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
3581a989b97SToby Isaac   for (i = 0, l = 0; i < n && l < k; i++) {
3591a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
3601a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
3611a989b97SToby Isaac 
3621a989b97SToby Isaac     if (j < Nminuskminus) {
3631a989b97SToby Isaac       subset[l++] = i;
3641a989b97SToby Isaac       Nk = Nminuskminus;
3651a989b97SToby Isaac     } else {
3661a989b97SToby Isaac       j -= Nminuskminus;
3671a989b97SToby Isaac       Nk = Nminusk;
3681a989b97SToby Isaac     }
3691a989b97SToby Isaac   }
3701a989b97SToby Isaac   PetscFunctionReturn(0);
3711a989b97SToby Isaac }
3721a989b97SToby Isaac 
373fad4db65SToby Isaac /*MC
374fad4db65SToby 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.
375fad4db65SToby Isaac 
376fad4db65SToby Isaac    Input Arguments:
377fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
378fad4db65SToby Isaac .  k - an integer in [0, n]
379fad4db65SToby Isaac -  subset - an ordered subset of the integers [0, ..., n - 1]
380fad4db65SToby Isaac 
381fad4db65SToby Isaac    Output Arguments:
382fad4db65SToby Isaac .  index - the rank of the subset in lexicographic order
383fad4db65SToby Isaac 
384fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
385fad4db65SToby Isaac 
386fad4db65SToby Isaac    Level: beginner
387fad4db65SToby Isaac 
388fad4db65SToby Isaac .seealso: PetscDTEnumSubset()
389fad4db65SToby Isaac M*/
3901a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
3911a989b97SToby Isaac {
3921a989b97SToby Isaac   PetscInt       i, j = 0, l, Nk;
3931a989b97SToby Isaac   PetscErrorCode ierr;
3941a989b97SToby Isaac 
39528222859SToby Isaac   PetscFunctionBegin;
39628222859SToby Isaac   *index = -1;
397fad4db65SToby Isaac   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
3981a989b97SToby Isaac   for (i = 0, l = 0; i < n && l < k; i++) {
3991a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
4001a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
4011a989b97SToby Isaac 
4021a989b97SToby Isaac     if (subset[l] == i) {
4031a989b97SToby Isaac       l++;
4041a989b97SToby Isaac       Nk = Nminuskminus;
4051a989b97SToby Isaac     } else {
4061a989b97SToby Isaac       j += Nminuskminus;
4071a989b97SToby Isaac       Nk = Nminusk;
4081a989b97SToby Isaac     }
4091a989b97SToby Isaac   }
4101a989b97SToby Isaac   *index = j;
4111a989b97SToby Isaac   PetscFunctionReturn(0);
4121a989b97SToby Isaac }
4131a989b97SToby Isaac 
414fad4db65SToby Isaac /*MC
41528222859SToby 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.
416fad4db65SToby Isaac 
417fad4db65SToby Isaac    Input Arguments:
418fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
419fad4db65SToby Isaac .  k - an integer in [0, n]
420fad4db65SToby Isaac -  j - an index in [0, n choose k)
421fad4db65SToby Isaac 
422fad4db65SToby Isaac    Output Arguments:
423fad4db65SToby Isaac +  perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set.
42428222859SToby Isaac -  isOdd - if not NULL, return whether perm is an even or odd permutation.
425fad4db65SToby Isaac 
426fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
427fad4db65SToby Isaac 
428fad4db65SToby Isaac    Level: beginner
429fad4db65SToby Isaac 
430fad4db65SToby Isaac .seealso: PetscDTEnumSubset(), PetscDTSubsetIndex()
431fad4db65SToby Isaac M*/
432fad4db65SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd)
4331a989b97SToby Isaac {
4341a989b97SToby Isaac   PetscInt       i, l, m, *subcomp, Nk;
4351a989b97SToby Isaac   PetscInt       odd;
4361a989b97SToby Isaac   PetscErrorCode ierr;
4371a989b97SToby Isaac 
43828222859SToby Isaac   PetscFunctionBegin;
43928222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
440fad4db65SToby Isaac   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
4411a989b97SToby Isaac   odd = 0;
442fad4db65SToby Isaac   subcomp = &perm[k];
4431a989b97SToby Isaac   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
4441a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
4451a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
4461a989b97SToby Isaac 
4471a989b97SToby Isaac     if (j < Nminuskminus) {
448fad4db65SToby Isaac       perm[l++] = i;
4491a989b97SToby Isaac       Nk = Nminuskminus;
4501a989b97SToby Isaac     } else {
4511a989b97SToby Isaac       subcomp[m++] = i;
4521a989b97SToby Isaac       j -= Nminuskminus;
4531a989b97SToby Isaac       odd ^= ((k - l) & 1);
4541a989b97SToby Isaac       Nk = Nminusk;
4551a989b97SToby Isaac     }
4561a989b97SToby Isaac   }
4571a989b97SToby Isaac   for (; i < n; i++) {
4581a989b97SToby Isaac     subcomp[m++] = i;
4591a989b97SToby Isaac   }
4601a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
4611a989b97SToby Isaac   PetscFunctionReturn(0);
4621a989b97SToby Isaac }
4631a989b97SToby Isaac 
464ef0bb6c7SMatthew G. Knepley struct _p_PetscTabulation {
465ef0bb6c7SMatthew G. Knepley   PetscInt    K;    /* Indicates a k-jet, namely tabulated derviatives up to order k */
466ef0bb6c7SMatthew G. Knepley   PetscInt    Nr;   /* THe number of tabulation replicas (often 1) */
467ef0bb6c7SMatthew G. Knepley   PetscInt    Np;   /* The number of tabulation points in a replica */
468ef0bb6c7SMatthew G. Knepley   PetscInt    Nb;   /* The number of functions tabulated */
469ef0bb6c7SMatthew G. Knepley   PetscInt    Nc;   /* The number of function components */
470ef0bb6c7SMatthew G. Knepley   PetscInt    cdim; /* The coordinate dimension */
471ef0bb6c7SMatthew G. Knepley   PetscReal **T;    /* The tabulation T[K] of functions and their derivatives
472ef0bb6c7SMatthew G. Knepley                        T[0] = B[Nr*Np][Nb][Nc]:             The basis function values at quadrature points
473ef0bb6c7SMatthew G. Knepley                        T[1] = D[Nr*Np][Nb][Nc][cdim]:       The basis function derivatives at quadrature points
474ef0bb6c7SMatthew G. Knepley                        T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */
475ef0bb6c7SMatthew G. Knepley };
476ef0bb6c7SMatthew G. Knepley typedef struct _p_PetscTabulation *PetscTabulation;
477ef0bb6c7SMatthew G. Knepley 
47837045ce4SJed Brown #endif
479