xref: /petsc/include/petscdt.h (revision 895212167bc0968de2d7bec19da4897a3fe677db)
137045ce4SJed Brown /*
237045ce4SJed Brown   Common tools for constructing discretizations
337045ce4SJed Brown */
4a4963045SJacob Faibussowitsch #pragma once
537045ce4SJed Brown 
637045ce4SJed Brown #include <petscsys.h>
74366bac7SMatthew G. Knepley #include <petscdmtypes.h>
84366bac7SMatthew G. Knepley #include <petscistypes.h>
937045ce4SJed Brown 
10ac09b921SBarry Smith /* SUBMANSEC = DT */
11ac09b921SBarry Smith 
122cd22861SMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCQUADRATURE_CLASSID;
132cd22861SMatthew G. Knepley 
1421454ff5SMatthew G. Knepley /*S
1516a05f60SBarry Smith   PetscQuadrature - Quadrature rule for numerical integration.
1621454ff5SMatthew G. Knepley 
17329bbf4eSMatthew G. Knepley   Level: beginner
1821454ff5SMatthew G. Knepley 
19db781477SPatrick Sanan .seealso: `PetscQuadratureCreate()`, `PetscQuadratureDestroy()`
2021454ff5SMatthew G. Knepley S*/
2121454ff5SMatthew G. Knepley typedef struct _p_PetscQuadrature *PetscQuadrature;
2221454ff5SMatthew G. Knepley 
238272889dSSatish Balay /*E
24916e780bShannah_mairs   PetscGaussLobattoLegendreCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights
258272889dSSatish Balay 
2616a05f60SBarry Smith   Values:
2716a05f60SBarry Smith +  `PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA` - compute the nodes via linear algebra
2816a05f60SBarry Smith -  `PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON`         - compute the nodes by solving a nonlinear equation with Newton's method
2916a05f60SBarry Smith 
308272889dSSatish Balay   Level: intermediate
318272889dSSatish Balay 
3216a05f60SBarry Smith .seealso: `PetscQuadrature`
338272889dSSatish Balay E*/
349371c9d4SSatish Balay typedef enum {
359371c9d4SSatish Balay   PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,
369371c9d4SSatish Balay   PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON
379371c9d4SSatish Balay } PetscGaussLobattoLegendreCreateType;
388272889dSSatish Balay 
39d4afb720SToby Isaac /*E
40d4afb720SToby Isaac   PetscDTNodeType - A description of strategies for generating nodes (both
41d4afb720SToby Isaac   quadrature nodes and nodes for Lagrange polynomials)
42d4afb720SToby Isaac 
4316a05f60SBarry Smith   Values:
4416a05f60SBarry Smith + `PETSCDTNODES_DEFAULT`     - Nodes chosen by PETSc
4516a05f60SBarry Smith . `PETSCDTNODES_GAUSSJACOBI` - Nodes at either Gauss-Jacobi or Gauss-Lobatto-Jacobi quadrature points
4616a05f60SBarry Smith . `PETSCDTNODES_EQUISPACED`  - Nodes equispaced either including the endpoints or excluding them
4716a05f60SBarry Smith - `PETSCDTNODES_TANHSINH`    - Nodes at Tanh-Sinh quadrature points
48d4afb720SToby Isaac 
4916a05f60SBarry Smith   Level: intermediate
50d4afb720SToby Isaac 
5187497f52SBarry Smith   Note:
5287497f52SBarry Smith   A `PetscDTNodeType` can be paired with a `PetscBool` to indicate whether
5387497f52SBarry Smith   the nodes include endpoints or not, and in the case of `PETSCDT_GAUSSJACOBI`
54d4afb720SToby Isaac   with exponents for the weight function.
55d4afb720SToby Isaac 
5616a05f60SBarry Smith .seealso: `PetscQuadrature`
57d4afb720SToby Isaac E*/
589371c9d4SSatish Balay typedef enum {
599371c9d4SSatish Balay   PETSCDTNODES_DEFAULT = -1,
609371c9d4SSatish Balay   PETSCDTNODES_GAUSSJACOBI,
619371c9d4SSatish Balay   PETSCDTNODES_EQUISPACED,
629371c9d4SSatish Balay   PETSCDTNODES_TANHSINH
639371c9d4SSatish Balay } PetscDTNodeType;
64d4afb720SToby Isaac 
65d3c69ad0SToby Isaac PETSC_EXTERN const char *const *const PetscDTNodeTypes;
66d3c69ad0SToby Isaac 
67d3c69ad0SToby Isaac /*E
68d3c69ad0SToby Isaac   PetscDTSimplexQuadratureType - A description of classes of quadrature rules for simplices
69d3c69ad0SToby Isaac 
7016a05f60SBarry Smith   Values:
7116a05f60SBarry Smith +  `PETSCDTSIMPLEXQUAD_DEFAULT` - Quadrature rule chosen by PETSc
7216a05f60SBarry Smith .  `PETSCDTSIMPLEXQUAD_CONIC`   - Quadrature rules constructed as
73d3c69ad0SToby Isaac                                   conically-warped tensor products of 1D
74d3c69ad0SToby Isaac                                   Gauss-Jacobi quadrature rules.  These are
75d3c69ad0SToby Isaac                                   explicitly computable in any dimension for any
76d3c69ad0SToby Isaac                                   degree, and the tensor-product structure can be
77d3c69ad0SToby Isaac                                   exploited by sum-factorization methods, but
78d3c69ad0SToby Isaac                                   they are not efficient in terms of nodes per
79d3c69ad0SToby Isaac                                   polynomial degree.
8016a05f60SBarry Smith -  `PETSCDTSIMPLEXQUAD_MINSYM`  - Quadrature rules that are fully symmetric
81d3c69ad0SToby Isaac                                   (symmetries of the simplex preserve the nodes
82d3c69ad0SToby Isaac                                   and weights) with minimal (or near minimal)
83d3c69ad0SToby Isaac                                   number of nodes.  In dimensions higher than 1
84d3c69ad0SToby Isaac                                   these are not simple to compute, so lookup
85d3c69ad0SToby Isaac                                   tables are used.
86d3c69ad0SToby Isaac 
8716a05f60SBarry Smith   Level: intermediate
8816a05f60SBarry Smith 
8916a05f60SBarry Smith .seealso: `PetscQuadrature`, `PetscDTSimplexQuadrature()`
90d3c69ad0SToby Isaac E*/
919371c9d4SSatish Balay typedef enum {
929371c9d4SSatish Balay   PETSCDTSIMPLEXQUAD_DEFAULT = -1,
939371c9d4SSatish Balay   PETSCDTSIMPLEXQUAD_CONIC   = 0,
949371c9d4SSatish Balay   PETSCDTSIMPLEXQUAD_MINSYM
959371c9d4SSatish Balay } PetscDTSimplexQuadratureType;
96d3c69ad0SToby Isaac 
97d3c69ad0SToby Isaac PETSC_EXTERN const char *const *const PetscDTSimplexQuadratureTypes;
98d4afb720SToby Isaac 
9921454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *);
100c9638911SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *);
1014366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetCellType(PetscQuadrature, DMPolytopeType *);
1024366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetCellType(PetscQuadrature, DMPolytopeType);
103bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt *);
104bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt);
105a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt *);
106a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt);
1074f9ab2b4SJed Brown PETSC_EXTERN PetscErrorCode PetscQuadratureEqual(PetscQuadrature, PetscQuadrature, PetscBool *);
108a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt *, PetscInt *, PetscInt *, const PetscReal *[], const PetscReal *[]);
109a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal[], const PetscReal[]);
11021454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer);
11121454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *);
112a0845e3aSMatthew G. Knepley 
1132df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTensorQuadratureCreate(PetscQuadrature, PetscQuadrature, PetscQuadrature *);
11489710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *);
1154366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureComputePermutations(PetscQuadrature, PetscInt *, IS *[]);
11689710940SMatthew G. Knepley 
117907761f8SToby Isaac PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *);
118907761f8SToby Isaac 
11937045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *);
120fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiNorm(PetscReal, PetscReal, PetscInt, PetscReal *);
12194e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt, PetscReal, PetscReal, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *);
122fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEvalJet(PetscReal, PetscReal, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]);
123fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPKDEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]);
124d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedSize(PetscInt, PetscInt, PetscInt, PetscInt *);
125d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscInt, PetscReal[]);
12637045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt, PetscReal, PetscReal, PetscReal *, PetscReal *);
12794e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *);
12894e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *);
129916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt, PetscGaussLobattoLegendreCreateType, PetscReal *, PetscReal *);
130194825f6SJed Brown PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
131a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
132e6a796c3SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
133d3c69ad0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTSimplexQuadrature(PetscInt, PetscInt, PetscDTSimplexQuadratureType, PetscQuadrature *);
1344366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTCreateDefaultQuadrature(DMPolytopeType, PetscInt, PetscQuadrature *, PetscQuadrature *);
13537045ce4SJed Brown 
136b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
137d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
138d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
139b3c0f97bSTom Klotz 
140916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *);
141916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
142916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
143916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
144916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
145916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
146916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
147916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
148916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
149916e780bShannah_mairs 
1502dce792eSToby Isaac /*MC
1512dce792eSToby Isaac   PETSC_FORM_DEGREE_UNDEFINED - Indicates that a field does not have
1522dce792eSToby Isaac   a well-defined form degree in exterior calculus.
1532dce792eSToby Isaac 
1542dce792eSToby Isaac   Level: advanced
1552dce792eSToby Isaac 
1562dce792eSToby Isaac .seealso: `PetscDTAltV`, `PetscDualSpaceGetFormDegree()`
1572dce792eSToby Isaac M*/
1582dce792eSToby Isaac #define PETSC_FORM_DEGREE_UNDEFINED PETSC_INT_MIN
1592dce792eSToby Isaac 
1601a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1611a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1621a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
1631a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
1641a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
1651a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1661a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
167dda711d0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
1681a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
1691a989b97SToby Isaac 
170d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt, PetscInt, const PetscInt[], PetscInt *);
171d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt, PetscInt, PetscInt, PetscInt[]);
172fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGradedOrderToIndex(PetscInt, const PetscInt[], PetscInt *);
173fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToGradedOrder(PetscInt, PetscInt, PetscInt[]);
174d4afb720SToby Isaac 
175fad4db65SToby Isaac #if defined(PETSC_USE_64BIT_INDICES)
176fad4db65SToby Isaac   #define PETSC_FACTORIAL_MAX 20
177fad4db65SToby Isaac   #define PETSC_BINOMIAL_MAX  61
178fad4db65SToby Isaac #else
179fad4db65SToby Isaac   #define PETSC_FACTORIAL_MAX 12
180fad4db65SToby Isaac   #define PETSC_BINOMIAL_MAX  29
181fad4db65SToby Isaac #endif
182fad4db65SToby Isaac 
183fad4db65SToby Isaac /*MC
184fad4db65SToby Isaac    PetscDTFactorial - Approximate n! as a real number
185fad4db65SToby Isaac 
1864165533cSJose E. Roman    Input Parameter:
187fad4db65SToby Isaac .  n - a non-negative integer
188fad4db65SToby Isaac 
1894165533cSJose E. Roman    Output Parameter:
190fad4db65SToby Isaac .  factorial - n!
191fad4db65SToby Isaac 
192fad4db65SToby Isaac    Level: beginner
19316a05f60SBarry Smith 
19416a05f60SBarry Smith .seealso: `PetscDTFactorialInt()`, `PetscDTBinomialInt()`, `PetscDTBinomial()`
195fad4db65SToby Isaac M*/
196d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial)
197d71ae5a4SJacob Faibussowitsch {
198fad4db65SToby Isaac   PetscReal f = 1.0;
199fad4db65SToby Isaac 
200fad4db65SToby Isaac   PetscFunctionBegin;
201e2ab39ccSLisandro Dalcin   *factorial = -1.0;
20263a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %" PetscInt_FMT, n);
2035f80ce2aSJacob Faibussowitsch   for (PetscInt i = 1; i < n + 1; ++i) f *= (PetscReal)i;
204fad4db65SToby Isaac   *factorial = f;
2053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
206fad4db65SToby Isaac }
207fad4db65SToby Isaac 
208fad4db65SToby Isaac /*MC
209fad4db65SToby Isaac    PetscDTFactorialInt - Compute n! as an integer
210fad4db65SToby Isaac 
2114165533cSJose E. Roman    Input Parameter:
212fad4db65SToby Isaac .  n - a non-negative integer
213fad4db65SToby Isaac 
2144165533cSJose E. Roman    Output Parameter:
215fad4db65SToby Isaac .  factorial - n!
216fad4db65SToby Isaac 
217fad4db65SToby Isaac    Level: beginner
218fad4db65SToby Isaac 
21916a05f60SBarry Smith    Note:
22095bd0b28SBarry Smith    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.
22116a05f60SBarry Smith 
22216a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTBinomialInt()`, `PetscDTBinomial()`
223fad4db65SToby Isaac M*/
224d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial)
225d71ae5a4SJacob Faibussowitsch {
226fad4db65SToby Isaac   PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600};
227fad4db65SToby Isaac 
22828222859SToby Isaac   PetscFunctionBegin;
22928222859SToby Isaac   *factorial = -1;
23063a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 0 && n <= PETSC_FACTORIAL_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of elements %" PetscInt_FMT " is not in supported range [0,%d]", n, PETSC_FACTORIAL_MAX);
231fad4db65SToby Isaac   if (n <= 12) {
232fad4db65SToby Isaac     *factorial = facLookup[n];
233fad4db65SToby Isaac   } else {
234fad4db65SToby Isaac     PetscInt f = facLookup[12];
235fad4db65SToby Isaac     PetscInt i;
236fad4db65SToby Isaac 
237fad4db65SToby Isaac     for (i = 13; i < n + 1; ++i) f *= i;
238fad4db65SToby Isaac     *factorial = f;
239fad4db65SToby Isaac   }
2403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
241fad4db65SToby Isaac }
242fad4db65SToby Isaac 
243fad4db65SToby Isaac /*MC
24495bd0b28SBarry Smith    PetscDTBinomial - Approximate the binomial coefficient `n` choose `k`
245fad4db65SToby Isaac 
2464165533cSJose E. Roman    Input Parameters:
247fad4db65SToby Isaac +  n - a non-negative integer
24895bd0b28SBarry Smith -  k - an integer between 0 and `n`, inclusive
249fad4db65SToby Isaac 
2504165533cSJose E. Roman    Output Parameter:
25195bd0b28SBarry Smith .  binomial - approximation of the binomial coefficient `n` choose `k`
252fad4db65SToby Isaac 
253fad4db65SToby Isaac    Level: beginner
25416a05f60SBarry Smith 
25516a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`
256fad4db65SToby Isaac M*/
257d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
258d71ae5a4SJacob Faibussowitsch {
2591a989b97SToby Isaac   PetscFunctionBeginHot;
260e2ab39ccSLisandro Dalcin   *binomial = -1.0;
26163a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 0 && k >= 0 && k <= n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial arguments (%" PetscInt_FMT " %" PetscInt_FMT ") must be non-negative, k <= n", n, k);
2621a989b97SToby Isaac   if (n <= 3) {
2639371c9d4SSatish Balay     PetscInt binomLookup[4][4] = {
2649371c9d4SSatish Balay       {1, 0, 0, 0},
2659371c9d4SSatish Balay       {1, 1, 0, 0},
2669371c9d4SSatish Balay       {1, 2, 1, 0},
2679371c9d4SSatish Balay       {1, 3, 3, 1}
2689371c9d4SSatish Balay     };
2691a989b97SToby Isaac 
270e2ab39ccSLisandro Dalcin     *binomial = (PetscReal)binomLookup[n][k];
2711a989b97SToby Isaac   } else {
272e2ab39ccSLisandro Dalcin     PetscReal binom = 1.0;
2731a989b97SToby Isaac 
2741a989b97SToby Isaac     k = PetscMin(k, n - k);
2755f80ce2aSJacob Faibussowitsch     for (PetscInt i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1);
2761a989b97SToby Isaac     *binomial = binom;
2771a989b97SToby Isaac   }
2783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2791a989b97SToby Isaac }
2801a989b97SToby Isaac 
281fad4db65SToby Isaac /*MC
28295bd0b28SBarry Smith    PetscDTBinomialInt - Compute the binomial coefficient `n` choose `k`
283fad4db65SToby Isaac 
28497bb3fdcSJose E. Roman    Input Parameters:
285fad4db65SToby Isaac +  n - a non-negative integer
28695bd0b28SBarry Smith -  k - an integer between 0 and `n`, inclusive
287fad4db65SToby Isaac 
28897bb3fdcSJose E. Roman    Output Parameter:
28995bd0b28SBarry Smith .  binomial - the binomial coefficient `n` choose `k`
290fad4db65SToby Isaac 
291fad4db65SToby Isaac    Level: beginner
29216a05f60SBarry Smith 
29316a05f60SBarry Smith    Note:
29416a05f60SBarry Smith    This is limited by integers that can be represented by `PetscInt`.
29516a05f60SBarry Smith 
29616a05f60SBarry Smith    Use `PetscDTBinomial()` for real number approximations of larger values
29716a05f60SBarry Smith 
29816a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTEnumPerm()`
299fad4db65SToby Isaac M*/
300d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial)
301d71ae5a4SJacob Faibussowitsch {
30228222859SToby Isaac   PetscInt bin;
30328222859SToby Isaac 
30428222859SToby Isaac   PetscFunctionBegin;
30528222859SToby Isaac   *binomial = -1;
30663a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 0 && k >= 0 && k <= n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial arguments (%" PetscInt_FMT " %" PetscInt_FMT ") must be non-negative, k <= n", n, k);
30763a3b9bcSJacob Faibussowitsch   PetscCheck(n <= PETSC_BINOMIAL_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial elements %" PetscInt_FMT " is larger than max for PetscInt, %d", n, PETSC_BINOMIAL_MAX);
308fad4db65SToby Isaac   if (n <= 3) {
3099371c9d4SSatish Balay     PetscInt binomLookup[4][4] = {
3109371c9d4SSatish Balay       {1, 0, 0, 0},
3119371c9d4SSatish Balay       {1, 1, 0, 0},
3129371c9d4SSatish Balay       {1, 2, 1, 0},
3139371c9d4SSatish Balay       {1, 3, 3, 1}
3149371c9d4SSatish Balay     };
315fad4db65SToby Isaac 
31628222859SToby Isaac     bin = binomLookup[n][k];
317fad4db65SToby Isaac   } else {
318fad4db65SToby Isaac     PetscInt binom = 1;
319fad4db65SToby Isaac 
320fad4db65SToby Isaac     k = PetscMin(k, n - k);
3215f80ce2aSJacob Faibussowitsch     for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
32228222859SToby Isaac     bin = binom;
323fad4db65SToby Isaac   }
32428222859SToby Isaac   *binomial = bin;
3253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
326fad4db65SToby Isaac }
327fad4db65SToby Isaac 
328fad4db65SToby Isaac /*MC
32916a05f60SBarry Smith    PetscDTEnumPerm - Get a permutation of `n` integers from its encoding into the integers [0, n!) as a sequence of swaps.
330fad4db65SToby Isaac 
3314165533cSJose E. Roman    Input Parameters:
332fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
3338cd1e013SToby Isaac -  k - an integer in [0, n!)
334fad4db65SToby Isaac 
3354165533cSJose E. Roman    Output Parameters:
336fad4db65SToby Isaac +  perm  - the permuted list of the integers [0, ..., n-1]
3372fe279fdSBarry Smith -  isOdd - if not `NULL`, returns whether the permutation used an even or odd number of swaps.
338fad4db65SToby Isaac 
33916a05f60SBarry Smith    Level: intermediate
340fad4db65SToby Isaac 
34116a05f60SBarry Smith    Notes:
34216a05f60SBarry Smith    A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation,
34316a05f60SBarry Smith    by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in
34416a05f60SBarry Smith    some position j >= i.  This swap is encoded as the difference (j - i).  The difference d_i at step i is less than
34516a05f60SBarry Smith    (n - i).  This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number
34616a05f60SBarry Smith    (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}.
34716a05f60SBarry Smith 
34816a05f60SBarry Smith    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.
34916a05f60SBarry Smith 
35016a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTPermIndex()`
351fad4db65SToby Isaac M*/
352d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd)
353d71ae5a4SJacob Faibussowitsch {
3541a989b97SToby Isaac   PetscInt  odd = 0;
3551a989b97SToby Isaac   PetscInt  i;
356fad4db65SToby Isaac   PetscInt  work[PETSC_FACTORIAL_MAX];
357fad4db65SToby Isaac   PetscInt *w;
3581a989b97SToby Isaac 
35928222859SToby Isaac   PetscFunctionBegin;
36028222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
36163a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 0 && n <= PETSC_FACTORIAL_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of elements %" PetscInt_FMT " is not in supported range [0,%d]", n, PETSC_FACTORIAL_MAX);
362810441c8SPierre Jolivet   if (n >= 2) {
363fad4db65SToby Isaac     w = &work[n - 2];
3641a989b97SToby Isaac     for (i = 2; i <= n; i++) {
3651a989b97SToby Isaac       *(w--) = k % i;
3661a989b97SToby Isaac       k /= i;
3671a989b97SToby Isaac     }
368810441c8SPierre Jolivet   }
3691a989b97SToby Isaac   for (i = 0; i < n; i++) perm[i] = i;
3701a989b97SToby Isaac   for (i = 0; i < n - 1; i++) {
3711a989b97SToby Isaac     PetscInt s    = work[i];
3721a989b97SToby Isaac     PetscInt swap = perm[i];
3731a989b97SToby Isaac 
3741a989b97SToby Isaac     perm[i]     = perm[i + s];
3751a989b97SToby Isaac     perm[i + s] = swap;
3761a989b97SToby Isaac     odd ^= (!!s);
3771a989b97SToby Isaac   }
3781a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
3793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3801a989b97SToby Isaac }
3811a989b97SToby Isaac 
382fad4db65SToby Isaac /*MC
38316a05f60SBarry Smith    PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!).  This inverts `PetscDTEnumPerm()`.
3848cd1e013SToby Isaac 
3854165533cSJose E. Roman    Input Parameters:
3868cd1e013SToby Isaac +  n    - a non-negative integer (see note about limits below)
3878cd1e013SToby Isaac -  perm - the permuted list of the integers [0, ..., n-1]
3888cd1e013SToby Isaac 
3894165533cSJose E. Roman    Output Parameters:
3908cd1e013SToby Isaac +  k     - an integer in [0, n!)
3912fe279fdSBarry Smith -  isOdd - if not `NULL`, returns whether the permutation used an even or odd number of swaps.
3928cd1e013SToby Isaac 
3938cd1e013SToby Isaac    Level: beginner
39416a05f60SBarry Smith 
39516a05f60SBarry Smith    Note:
39616a05f60SBarry Smith    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.
39716a05f60SBarry Smith 
39816a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`
3998cd1e013SToby Isaac M*/
400d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd)
401d71ae5a4SJacob Faibussowitsch {
4028cd1e013SToby Isaac   PetscInt odd = 0;
4038cd1e013SToby Isaac   PetscInt i, idx;
4048cd1e013SToby Isaac   PetscInt work[PETSC_FACTORIAL_MAX];
4058cd1e013SToby Isaac   PetscInt iwork[PETSC_FACTORIAL_MAX];
4068cd1e013SToby Isaac 
4078cd1e013SToby Isaac   PetscFunctionBeginHot;
40828222859SToby Isaac   *k = -1;
40928222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
41063a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 0 && n <= PETSC_FACTORIAL_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of elements %" PetscInt_FMT " is not in supported range [0,%d]", n, PETSC_FACTORIAL_MAX);
4118cd1e013SToby Isaac   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
4128cd1e013SToby Isaac   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
4138cd1e013SToby Isaac   for (idx = 0, i = 0; i < n - 1; i++) {
4148cd1e013SToby Isaac     PetscInt j    = perm[i];
4158cd1e013SToby Isaac     PetscInt icur = work[i];
4168cd1e013SToby Isaac     PetscInt jloc = iwork[j];
4178cd1e013SToby Isaac     PetscInt diff = jloc - i;
4188cd1e013SToby Isaac 
4198cd1e013SToby Isaac     idx = idx * (n - i) + diff;
4208cd1e013SToby Isaac     /* swap (i, jloc) */
4218cd1e013SToby Isaac     work[i]     = j;
4228cd1e013SToby Isaac     work[jloc]  = icur;
4238cd1e013SToby Isaac     iwork[j]    = i;
4248cd1e013SToby Isaac     iwork[icur] = jloc;
4258cd1e013SToby Isaac     odd ^= (!!diff);
4268cd1e013SToby Isaac   }
4278cd1e013SToby Isaac   *k = idx;
4288cd1e013SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
4293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4308cd1e013SToby Isaac }
4318cd1e013SToby Isaac 
4328cd1e013SToby Isaac /*MC
433fad4db65SToby Isaac    PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k).
434fad4db65SToby Isaac    The encoding is in lexicographic order.
435fad4db65SToby Isaac 
4364165533cSJose E. Roman    Input Parameters:
437fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
438fad4db65SToby Isaac .  k - an integer in [0, n]
439fad4db65SToby Isaac -  j - an index in [0, n choose k)
440fad4db65SToby Isaac 
4414165533cSJose E. Roman    Output Parameter:
442fad4db65SToby Isaac .  subset - the jth subset of size k of the integers [0, ..., n - 1]
443fad4db65SToby Isaac 
444fad4db65SToby Isaac    Level: beginner
445fad4db65SToby Isaac 
44616a05f60SBarry Smith    Note:
44716a05f60SBarry Smith    Limited by arguments such that `n` choose `k` can be represented by `PetscInt`
44816a05f60SBarry Smith 
44916a05f60SBarry Smith .seealso: `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()`
450fad4db65SToby Isaac M*/
451d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
452d71ae5a4SJacob Faibussowitsch {
4535f80ce2aSJacob Faibussowitsch   PetscInt Nk;
4541a989b97SToby Isaac 
4551a989b97SToby Isaac   PetscFunctionBeginHot;
4569566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
4575f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
4581a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
4591a989b97SToby Isaac     PetscInt Nminusk      = Nk - Nminuskminus;
4601a989b97SToby Isaac 
4611a989b97SToby Isaac     if (j < Nminuskminus) {
4621a989b97SToby Isaac       subset[l++] = i;
4631a989b97SToby Isaac       Nk          = Nminuskminus;
4641a989b97SToby Isaac     } else {
4651a989b97SToby Isaac       j -= Nminuskminus;
4661a989b97SToby Isaac       Nk = Nminusk;
4671a989b97SToby Isaac     }
4681a989b97SToby Isaac   }
4693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4701a989b97SToby Isaac }
4711a989b97SToby Isaac 
472fad4db65SToby Isaac /*MC
47387497f52SBarry Smith    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.
47487497f52SBarry Smith    This is the inverse of `PetscDTEnumSubset`.
475fad4db65SToby Isaac 
4764165533cSJose E. Roman    Input Parameters:
477fad4db65SToby Isaac +  n      - a non-negative integer (see note about limits below)
478fad4db65SToby Isaac .  k      - an integer in [0, n]
479fad4db65SToby Isaac -  subset - an ordered subset of the integers [0, ..., n - 1]
480fad4db65SToby Isaac 
4814165533cSJose E. Roman    Output Parameter:
482fad4db65SToby Isaac .  index - the rank of the subset in lexicographic order
483fad4db65SToby Isaac 
484fad4db65SToby Isaac    Level: beginner
485fad4db65SToby Isaac 
48616a05f60SBarry Smith    Note:
48716a05f60SBarry Smith    Limited by arguments such that `n` choose `k` can be represented by `PetscInt`
48816a05f60SBarry Smith 
48916a05f60SBarry Smith .seealso: `PetscDTEnumSubset()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()`
490fad4db65SToby Isaac M*/
491d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
492d71ae5a4SJacob Faibussowitsch {
4935f80ce2aSJacob Faibussowitsch   PetscInt j = 0, Nk;
4941a989b97SToby Isaac 
49528222859SToby Isaac   PetscFunctionBegin;
49628222859SToby Isaac   *index = -1;
4979566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
4985f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
4991a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
5001a989b97SToby Isaac     PetscInt Nminusk      = Nk - Nminuskminus;
5011a989b97SToby Isaac 
5021a989b97SToby Isaac     if (subset[l] == i) {
5031a989b97SToby Isaac       l++;
5041a989b97SToby Isaac       Nk = Nminuskminus;
5051a989b97SToby Isaac     } else {
5061a989b97SToby Isaac       j += Nminuskminus;
5071a989b97SToby Isaac       Nk = Nminusk;
5081a989b97SToby Isaac     }
5091a989b97SToby Isaac   }
5101a989b97SToby Isaac   *index = j;
5113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5121a989b97SToby Isaac }
5131a989b97SToby Isaac 
514fad4db65SToby Isaac /*MC
515*89521216SMatthew G. Knepley    PetscDTEnumSplit - 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.
516fad4db65SToby Isaac 
5174165533cSJose E. Roman    Input Parameters:
518fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
519fad4db65SToby Isaac .  k - an integer in [0, n]
520fad4db65SToby Isaac -  j - an index in [0, n choose k)
521fad4db65SToby Isaac 
5224165533cSJose E. Roman    Output Parameters:
523fad4db65SToby Isaac +  perm  - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set.
5242fe279fdSBarry Smith -  isOdd - if not `NULL`, return whether perm is an even or odd permutation.
525fad4db65SToby Isaac 
526fad4db65SToby Isaac    Level: beginner
527fad4db65SToby Isaac 
52816a05f60SBarry Smith    Note:
52916a05f60SBarry Smith    Limited by arguments such that `n` choose `k` can be represented by `PetscInt`
53016a05f60SBarry Smith 
53116a05f60SBarry Smith .seealso: `PetscDTEnumSubset()`, `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`,
53216a05f60SBarry Smith           `PetscDTPermIndex()`
533fad4db65SToby Isaac M*/
534d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd)
535d71ae5a4SJacob Faibussowitsch {
5365f80ce2aSJacob Faibussowitsch   PetscInt  i, l, m, Nk, odd = 0;
5378e3a54c0SPierre Jolivet   PetscInt *subcomp = PetscSafePointerPlusOffset(perm, k);
5381a989b97SToby Isaac 
53928222859SToby Isaac   PetscFunctionBegin;
54028222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
5419566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
5421a989b97SToby Isaac   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
5431a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
5441a989b97SToby Isaac     PetscInt Nminusk      = Nk - Nminuskminus;
5451a989b97SToby Isaac 
5461a989b97SToby Isaac     if (j < Nminuskminus) {
547fad4db65SToby Isaac       perm[l++] = i;
5481a989b97SToby Isaac       Nk        = Nminuskminus;
5491a989b97SToby Isaac     } else {
5501a989b97SToby Isaac       subcomp[m++] = i;
5511a989b97SToby Isaac       j -= Nminuskminus;
5521a989b97SToby Isaac       odd ^= ((k - l) & 1);
5531a989b97SToby Isaac       Nk = Nminusk;
5541a989b97SToby Isaac     }
5551a989b97SToby Isaac   }
5565f80ce2aSJacob Faibussowitsch   for (; i < n; i++) subcomp[m++] = i;
5571a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
5583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5591a989b97SToby Isaac }
5601a989b97SToby Isaac 
561ef0bb6c7SMatthew G. Knepley struct _p_PetscTabulation {
562a5b23f4aSJose E. Roman   PetscInt    K;    /* Indicates a k-jet, namely tabulated derivatives up to order k */
56319815104SMartin Diehl   PetscInt    Nr;   /* The number of tabulation replicas (often 1) */
564ef0bb6c7SMatthew G. Knepley   PetscInt    Np;   /* The number of tabulation points in a replica */
565ef0bb6c7SMatthew G. Knepley   PetscInt    Nb;   /* The number of functions tabulated */
566ef0bb6c7SMatthew G. Knepley   PetscInt    Nc;   /* The number of function components */
567ef0bb6c7SMatthew G. Knepley   PetscInt    cdim; /* The coordinate dimension */
568ef0bb6c7SMatthew G. Knepley   PetscReal **T;    /* The tabulation T[K] of functions and their derivatives
569ef0bb6c7SMatthew G. Knepley                        T[0] = B[Nr*Np][Nb][Nc]:             The basis function values at quadrature points
570ef0bb6c7SMatthew G. Knepley                        T[1] = D[Nr*Np][Nb][Nc][cdim]:       The basis function derivatives at quadrature points
571ef0bb6c7SMatthew G. Knepley                        T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */
572ef0bb6c7SMatthew G. Knepley };
573ef0bb6c7SMatthew G. Knepley typedef struct _p_PetscTabulation *PetscTabulation;
574ef0bb6c7SMatthew G. Knepley 
575d6685f55SMatthew G. Knepley typedef PetscErrorCode (*PetscProbFunc)(const PetscReal[], const PetscReal[], PetscReal[]);
576d6685f55SMatthew G. Knepley 
5779371c9d4SSatish Balay typedef enum {
5789371c9d4SSatish Balay   DTPROB_DENSITY_CONSTANT,
5799371c9d4SSatish Balay   DTPROB_DENSITY_GAUSSIAN,
5809371c9d4SSatish Balay   DTPROB_DENSITY_MAXWELL_BOLTZMANN,
5819371c9d4SSatish Balay   DTPROB_NUM_DENSITY
5829371c9d4SSatish Balay } DTProbDensityType;
583d6685f55SMatthew G. Knepley PETSC_EXTERN const char *const DTProbDensityTypes[];
584d6685f55SMatthew G. Knepley 
585d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
586d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
587d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
588d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
589d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
590d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
591d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
592d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
593d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
594d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
595d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
596ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]);
597ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]);
598d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
599d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
600d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
601ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
602ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
603ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
604ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
605ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
606ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
607d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFunc *, PetscProbFunc *, PetscProbFunc *);
608d6685f55SMatthew G. Knepley 
609d6685f55SMatthew G. Knepley #include <petscvec.h>
610d6685f55SMatthew G. Knepley 
611d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFunc, PetscReal *);
612