xref: /petsc/include/petscdt.h (revision 4366bac74fcc7971b912b77ced138eb774c2c9eb)
137045ce4SJed Brown /*
237045ce4SJed Brown   Common tools for constructing discretizations
337045ce4SJed Brown */
46524c165SJacob Faibussowitsch #ifndef PETSCDT_H
526bd1501SBarry Smith #define PETSCDT_H
637045ce4SJed Brown 
737045ce4SJed Brown #include <petscsys.h>
8*4366bac7SMatthew G. Knepley #include <petscdmtypes.h>
9*4366bac7SMatthew G. Knepley #include <petscistypes.h>
1037045ce4SJed Brown 
11ac09b921SBarry Smith /* SUBMANSEC = DT */
12ac09b921SBarry Smith 
132cd22861SMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCQUADRATURE_CLASSID;
142cd22861SMatthew G. Knepley 
1521454ff5SMatthew G. Knepley /*S
1621454ff5SMatthew G. Knepley   PetscQuadrature - Quadrature rule for integration.
1721454ff5SMatthew G. Knepley 
18329bbf4eSMatthew G. Knepley   Level: beginner
1921454ff5SMatthew G. Knepley 
20db781477SPatrick Sanan .seealso: `PetscQuadratureCreate()`, `PetscQuadratureDestroy()`
2121454ff5SMatthew G. Knepley S*/
2221454ff5SMatthew G. Knepley typedef struct _p_PetscQuadrature *PetscQuadrature;
2321454ff5SMatthew G. Knepley 
248272889dSSatish Balay /*E
25916e780bShannah_mairs   PetscGaussLobattoLegendreCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights
268272889dSSatish Balay 
278272889dSSatish Balay   Level: intermediate
288272889dSSatish Balay 
2987497f52SBarry Smith $  `PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA` - compute the nodes via linear algebra
3087497f52SBarry Smith $  `PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON` - compute the nodes by solving a nonlinear equation with Newton's method
318272889dSSatish Balay 
328272889dSSatish Balay E*/
339371c9d4SSatish Balay typedef enum {
349371c9d4SSatish Balay   PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,
359371c9d4SSatish Balay   PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON
369371c9d4SSatish Balay } PetscGaussLobattoLegendreCreateType;
378272889dSSatish Balay 
38d4afb720SToby Isaac /*E
39d4afb720SToby Isaac   PetscDTNodeType - A description of strategies for generating nodes (both
40d4afb720SToby Isaac   quadrature nodes and nodes for Lagrange polynomials)
41d4afb720SToby Isaac 
42d4afb720SToby Isaac   Level: intermediate
43d4afb720SToby Isaac 
4487497f52SBarry Smith $  `PETSCDTNODES_DEFAULT` - Nodes chosen by PETSc
4587497f52SBarry Smith $  `PETSCDTNODES_GAUSSJACOBI` - Nodes at either Gauss-Jacobi or Gauss-Lobatto-Jacobi quadrature points
4687497f52SBarry Smith $  `PETSCDTNODES_EQUISPACED` - Nodes equispaced either including the endpoints or excluding them
4787497f52SBarry Smith $  `PETSCDTNODES_TANHSINH` - Nodes at Tanh-Sinh quadrature points
48d4afb720SToby Isaac 
4987497f52SBarry Smith   Note:
5087497f52SBarry Smith   A `PetscDTNodeType` can be paired with a `PetscBool` to indicate whether
5187497f52SBarry Smith   the nodes include endpoints or not, and in the case of `PETSCDT_GAUSSJACOBI`
52d4afb720SToby Isaac   with exponents for the weight function.
53d4afb720SToby Isaac 
54d4afb720SToby Isaac E*/
559371c9d4SSatish Balay typedef enum {
569371c9d4SSatish Balay   PETSCDTNODES_DEFAULT = -1,
579371c9d4SSatish Balay   PETSCDTNODES_GAUSSJACOBI,
589371c9d4SSatish Balay   PETSCDTNODES_EQUISPACED,
599371c9d4SSatish Balay   PETSCDTNODES_TANHSINH
609371c9d4SSatish Balay } PetscDTNodeType;
61d4afb720SToby Isaac 
62d3c69ad0SToby Isaac PETSC_EXTERN const char *const *const PetscDTNodeTypes;
63d3c69ad0SToby Isaac 
64d3c69ad0SToby Isaac /*E
65d3c69ad0SToby Isaac   PetscDTSimplexQuadratureType - A description of classes of quadrature rules for simplices
66d3c69ad0SToby Isaac 
67d3c69ad0SToby Isaac   Level: intermediate
68d3c69ad0SToby Isaac 
6987497f52SBarry Smith $  `PETSCDTSIMPLEXQUAD_DEFAULT` - Quadrature rule chosen by PETSc
7087497f52SBarry Smith $  `PETSCDTSIMPLEXQUAD_CONIC`   - Quadrature rules constructed as
71d3c69ad0SToby Isaac                                 conically-warped tensor products of 1D
72d3c69ad0SToby Isaac                                 Gauss-Jacobi quadrature rules.  These are
73d3c69ad0SToby Isaac                                 explicitly computable in any dimension for any
74d3c69ad0SToby Isaac                                 degree, and the tensor-product structure can be
75d3c69ad0SToby Isaac                                 exploited by sum-factorization methods, but
76d3c69ad0SToby Isaac                                 they are not efficient in terms of nodes per
77d3c69ad0SToby Isaac                                 polynomial degree.
7887497f52SBarry Smith $  `PETSCDTSIMPLEXQUAD_MINSYM`  - Quadrature rules that are fully symmetric
79d3c69ad0SToby Isaac                                 (symmetries of the simplex preserve the nodes
80d3c69ad0SToby Isaac                                 and weights) with minimal (or near minimal)
81d3c69ad0SToby Isaac                                 number of nodes.  In dimensions higher than 1
82d3c69ad0SToby Isaac                                 these are not simple to compute, so lookup
83d3c69ad0SToby Isaac                                 tables are used.
84d3c69ad0SToby Isaac 
85d3c69ad0SToby Isaac .seealso: `PetscDTSimplexQuadrature()`
86d3c69ad0SToby Isaac E*/
879371c9d4SSatish Balay typedef enum {
889371c9d4SSatish Balay   PETSCDTSIMPLEXQUAD_DEFAULT = -1,
899371c9d4SSatish Balay   PETSCDTSIMPLEXQUAD_CONIC   = 0,
909371c9d4SSatish Balay   PETSCDTSIMPLEXQUAD_MINSYM
919371c9d4SSatish Balay } PetscDTSimplexQuadratureType;
92d3c69ad0SToby Isaac 
93d3c69ad0SToby Isaac PETSC_EXTERN const char *const *const PetscDTSimplexQuadratureTypes;
94d4afb720SToby Isaac 
9521454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *);
96c9638911SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *);
97*4366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetCellType(PetscQuadrature, DMPolytopeType *);
98*4366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetCellType(PetscQuadrature, DMPolytopeType);
99bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt *);
100bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt);
101a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt *);
102a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt);
1034f9ab2b4SJed Brown PETSC_EXTERN PetscErrorCode PetscQuadratureEqual(PetscQuadrature, PetscQuadrature, PetscBool *);
104a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt *, PetscInt *, PetscInt *, const PetscReal *[], const PetscReal *[]);
105a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal[], const PetscReal[]);
10621454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer);
10721454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *);
108a0845e3aSMatthew G. Knepley 
1092df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTensorQuadratureCreate(PetscQuadrature, PetscQuadrature, PetscQuadrature *);
11089710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *);
111*4366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureComputePermutations(PetscQuadrature, PetscInt *, IS *[]);
11289710940SMatthew G. Knepley 
113907761f8SToby Isaac PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *);
114907761f8SToby Isaac 
11537045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *);
116fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiNorm(PetscReal, PetscReal, PetscInt, PetscReal *);
11794e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt, PetscReal, PetscReal, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *);
118fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEvalJet(PetscReal, PetscReal, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]);
119fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPKDEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]);
120d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedSize(PetscInt, PetscInt, PetscInt, PetscInt *);
121d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscInt, PetscReal[]);
12237045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt, PetscReal, PetscReal, PetscReal *, PetscReal *);
12394e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *);
12494e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *);
125916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt, PetscGaussLobattoLegendreCreateType, PetscReal *, PetscReal *);
126194825f6SJed Brown PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
127a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
128e6a796c3SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
129d3c69ad0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTSimplexQuadrature(PetscInt, PetscInt, PetscDTSimplexQuadratureType, PetscQuadrature *);
130*4366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTCreateDefaultQuadrature(DMPolytopeType, PetscInt, PetscQuadrature *, PetscQuadrature *);
13137045ce4SJed Brown 
132b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
133d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
134d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
135b3c0f97bSTom Klotz 
136916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *);
137916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
138916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
139916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
140916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
141916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
142916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
143916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
144916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
145916e780bShannah_mairs 
1461a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1471a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1481a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
1491a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
1501a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
1511a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1521a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
153dda711d0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
1541a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
1551a989b97SToby Isaac 
156d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt, PetscInt, const PetscInt[], PetscInt *);
157d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt, PetscInt, PetscInt, PetscInt[]);
158fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGradedOrderToIndex(PetscInt, const PetscInt[], PetscInt *);
159fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToGradedOrder(PetscInt, PetscInt, PetscInt[]);
160d4afb720SToby Isaac 
161fad4db65SToby Isaac #if defined(PETSC_USE_64BIT_INDICES)
162fad4db65SToby Isaac   #define PETSC_FACTORIAL_MAX 20
163fad4db65SToby Isaac   #define PETSC_BINOMIAL_MAX  61
164fad4db65SToby Isaac #else
165fad4db65SToby Isaac   #define PETSC_FACTORIAL_MAX 12
166fad4db65SToby Isaac   #define PETSC_BINOMIAL_MAX  29
167fad4db65SToby Isaac #endif
168fad4db65SToby Isaac 
169fad4db65SToby Isaac /*MC
170fad4db65SToby Isaac    PetscDTFactorial - Approximate n! as a real number
171fad4db65SToby Isaac 
1724165533cSJose E. Roman    Input Parameter:
173fad4db65SToby Isaac .  n - a non-negative integer
174fad4db65SToby Isaac 
1754165533cSJose E. Roman    Output Parameter:
176fad4db65SToby Isaac .  factorial - n!
177fad4db65SToby Isaac 
178fad4db65SToby Isaac    Level: beginner
179fad4db65SToby Isaac M*/
180d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial)
181d71ae5a4SJacob Faibussowitsch {
182fad4db65SToby Isaac   PetscReal f = 1.0;
183fad4db65SToby Isaac 
184fad4db65SToby Isaac   PetscFunctionBegin;
185e2ab39ccSLisandro Dalcin   *factorial = -1.0;
18663a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %" PetscInt_FMT, n);
1875f80ce2aSJacob Faibussowitsch   for (PetscInt i = 1; i < n + 1; ++i) f *= (PetscReal)i;
188fad4db65SToby Isaac   *factorial = f;
1893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
190fad4db65SToby Isaac }
191fad4db65SToby Isaac 
192fad4db65SToby Isaac /*MC
193fad4db65SToby Isaac    PetscDTFactorialInt - Compute n! as an integer
194fad4db65SToby Isaac 
1954165533cSJose E. Roman    Input Parameter:
196fad4db65SToby Isaac .  n - a non-negative integer
197fad4db65SToby Isaac 
1984165533cSJose E. Roman    Output Parameter:
199fad4db65SToby Isaac .  factorial - n!
200fad4db65SToby Isaac 
201fad4db65SToby Isaac    Level: beginner
202fad4db65SToby Isaac 
203fad4db65SToby 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.
204fad4db65SToby Isaac M*/
205d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial)
206d71ae5a4SJacob Faibussowitsch {
207fad4db65SToby Isaac   PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600};
208fad4db65SToby Isaac 
20928222859SToby Isaac   PetscFunctionBegin;
21028222859SToby Isaac   *factorial = -1;
21163a3b9bcSJacob 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);
212fad4db65SToby Isaac   if (n <= 12) {
213fad4db65SToby Isaac     *factorial = facLookup[n];
214fad4db65SToby Isaac   } else {
215fad4db65SToby Isaac     PetscInt f = facLookup[12];
216fad4db65SToby Isaac     PetscInt i;
217fad4db65SToby Isaac 
218fad4db65SToby Isaac     for (i = 13; i < n + 1; ++i) f *= i;
219fad4db65SToby Isaac     *factorial = f;
220fad4db65SToby Isaac   }
2213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
222fad4db65SToby Isaac }
223fad4db65SToby Isaac 
224fad4db65SToby Isaac /*MC
225fad4db65SToby Isaac    PetscDTBinomial - Approximate the binomial coefficient "n choose k"
226fad4db65SToby Isaac 
2274165533cSJose E. Roman    Input Parameters:
228fad4db65SToby Isaac +  n - a non-negative integer
229fad4db65SToby Isaac -  k - an integer between 0 and n, inclusive
230fad4db65SToby Isaac 
2314165533cSJose E. Roman    Output Parameter:
232fad4db65SToby Isaac .  binomial - approximation of the binomial coefficient n choose k
233fad4db65SToby Isaac 
234fad4db65SToby Isaac    Level: beginner
235fad4db65SToby Isaac M*/
236d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
237d71ae5a4SJacob Faibussowitsch {
2381a989b97SToby Isaac   PetscFunctionBeginHot;
239e2ab39ccSLisandro Dalcin   *binomial = -1.0;
24063a3b9bcSJacob 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);
2411a989b97SToby Isaac   if (n <= 3) {
2429371c9d4SSatish Balay     PetscInt binomLookup[4][4] = {
2439371c9d4SSatish Balay       {1, 0, 0, 0},
2449371c9d4SSatish Balay       {1, 1, 0, 0},
2459371c9d4SSatish Balay       {1, 2, 1, 0},
2469371c9d4SSatish Balay       {1, 3, 3, 1}
2479371c9d4SSatish Balay     };
2481a989b97SToby Isaac 
249e2ab39ccSLisandro Dalcin     *binomial = (PetscReal)binomLookup[n][k];
2501a989b97SToby Isaac   } else {
251e2ab39ccSLisandro Dalcin     PetscReal binom = 1.0;
2521a989b97SToby Isaac 
2531a989b97SToby Isaac     k = PetscMin(k, n - k);
2545f80ce2aSJacob Faibussowitsch     for (PetscInt i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1);
2551a989b97SToby Isaac     *binomial = binom;
2561a989b97SToby Isaac   }
2573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2581a989b97SToby Isaac }
2591a989b97SToby Isaac 
260fad4db65SToby Isaac /*MC
261fad4db65SToby Isaac    PetscDTBinomialInt - Compute the binomial coefficient "n choose k"
262fad4db65SToby Isaac 
26397bb3fdcSJose E. Roman    Input Parameters:
264fad4db65SToby Isaac +  n - a non-negative integer
265fad4db65SToby Isaac -  k - an integer between 0 and n, inclusive
266fad4db65SToby Isaac 
26797bb3fdcSJose E. Roman    Output Parameter:
268fad4db65SToby Isaac .  binomial - the binomial coefficient n choose k
269fad4db65SToby Isaac 
27087497f52SBarry Smith    Note: this is limited by integers that can be represented by `PetscInt`
271fad4db65SToby Isaac 
272fad4db65SToby Isaac    Level: beginner
273fad4db65SToby Isaac M*/
274d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial)
275d71ae5a4SJacob Faibussowitsch {
27628222859SToby Isaac   PetscInt bin;
27728222859SToby Isaac 
27828222859SToby Isaac   PetscFunctionBegin;
27928222859SToby Isaac   *binomial = -1;
28063a3b9bcSJacob 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);
28163a3b9bcSJacob 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);
282fad4db65SToby Isaac   if (n <= 3) {
2839371c9d4SSatish Balay     PetscInt binomLookup[4][4] = {
2849371c9d4SSatish Balay       {1, 0, 0, 0},
2859371c9d4SSatish Balay       {1, 1, 0, 0},
2869371c9d4SSatish Balay       {1, 2, 1, 0},
2879371c9d4SSatish Balay       {1, 3, 3, 1}
2889371c9d4SSatish Balay     };
289fad4db65SToby Isaac 
29028222859SToby Isaac     bin = binomLookup[n][k];
291fad4db65SToby Isaac   } else {
292fad4db65SToby Isaac     PetscInt binom = 1;
293fad4db65SToby Isaac 
294fad4db65SToby Isaac     k = PetscMin(k, n - k);
2955f80ce2aSJacob Faibussowitsch     for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
29628222859SToby Isaac     bin = binom;
297fad4db65SToby Isaac   }
29828222859SToby Isaac   *binomial = bin;
2993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
300fad4db65SToby Isaac }
301fad4db65SToby Isaac 
302fad4db65SToby Isaac /*MC
303fad4db65SToby Isaac    PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps.
304fad4db65SToby Isaac 
305fad4db65SToby Isaac    A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation,
306fad4db65SToby Isaac    by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in
30728222859SToby Isaac    some position j >= i.  This swap is encoded as the difference (j - i).  The difference d_i at step i is less than
30828222859SToby Isaac    (n - i).  This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number
309fad4db65SToby Isaac    (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}.
310fad4db65SToby Isaac 
3114165533cSJose E. Roman    Input Parameters:
312fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
3138cd1e013SToby Isaac -  k - an integer in [0, n!)
314fad4db65SToby Isaac 
3154165533cSJose E. Roman    Output Parameters:
316fad4db65SToby Isaac +  perm - the permuted list of the integers [0, ..., n-1]
317d5b43468SJose E. Roman -  isOdd - if not NULL, returns whether the permutation used an even or odd number of swaps.
318fad4db65SToby Isaac 
31987497f52SBarry Smith    Note:
32087497f52SBarry 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.
321fad4db65SToby Isaac 
322fad4db65SToby Isaac    Level: beginner
323fad4db65SToby Isaac M*/
324d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd)
325d71ae5a4SJacob Faibussowitsch {
3261a989b97SToby Isaac   PetscInt  odd = 0;
3271a989b97SToby Isaac   PetscInt  i;
328fad4db65SToby Isaac   PetscInt  work[PETSC_FACTORIAL_MAX];
329fad4db65SToby Isaac   PetscInt *w;
3301a989b97SToby Isaac 
33128222859SToby Isaac   PetscFunctionBegin;
33228222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
33363a3b9bcSJacob 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);
334fad4db65SToby Isaac   w = &work[n - 2];
3351a989b97SToby Isaac   for (i = 2; i <= n; i++) {
3361a989b97SToby Isaac     *(w--) = k % i;
3371a989b97SToby Isaac     k /= i;
3381a989b97SToby Isaac   }
3391a989b97SToby Isaac   for (i = 0; i < n; i++) perm[i] = i;
3401a989b97SToby Isaac   for (i = 0; i < n - 1; i++) {
3411a989b97SToby Isaac     PetscInt s    = work[i];
3421a989b97SToby Isaac     PetscInt swap = perm[i];
3431a989b97SToby Isaac 
3441a989b97SToby Isaac     perm[i]     = perm[i + s];
3451a989b97SToby Isaac     perm[i + s] = swap;
3461a989b97SToby Isaac     odd ^= (!!s);
3471a989b97SToby Isaac   }
3481a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
3493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3501a989b97SToby Isaac }
3511a989b97SToby Isaac 
352fad4db65SToby Isaac /*MC
35387497f52SBarry Smith    PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!).  This inverts `PetscDTEnumPerm`.
3548cd1e013SToby Isaac 
3554165533cSJose E. Roman    Input Parameters:
3568cd1e013SToby Isaac +  n - a non-negative integer (see note about limits below)
3578cd1e013SToby Isaac -  perm - the permuted list of the integers [0, ..., n-1]
3588cd1e013SToby Isaac 
3594165533cSJose E. Roman    Output Parameters:
3608cd1e013SToby Isaac +  k - an integer in [0, n!)
361d5b43468SJose E. Roman -  isOdd - if not NULL, returns whether the permutation used an even or odd number of swaps.
3628cd1e013SToby Isaac 
36387497f52SBarry Smith    Note:
36487497f52SBarry 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.
3658cd1e013SToby Isaac 
3668cd1e013SToby Isaac    Level: beginner
3678cd1e013SToby Isaac M*/
368d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd)
369d71ae5a4SJacob Faibussowitsch {
3708cd1e013SToby Isaac   PetscInt odd = 0;
3718cd1e013SToby Isaac   PetscInt i, idx;
3728cd1e013SToby Isaac   PetscInt work[PETSC_FACTORIAL_MAX];
3738cd1e013SToby Isaac   PetscInt iwork[PETSC_FACTORIAL_MAX];
3748cd1e013SToby Isaac 
3758cd1e013SToby Isaac   PetscFunctionBeginHot;
37628222859SToby Isaac   *k = -1;
37728222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
37863a3b9bcSJacob 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);
3798cd1e013SToby Isaac   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
3808cd1e013SToby Isaac   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
3818cd1e013SToby Isaac   for (idx = 0, i = 0; i < n - 1; i++) {
3828cd1e013SToby Isaac     PetscInt j    = perm[i];
3838cd1e013SToby Isaac     PetscInt icur = work[i];
3848cd1e013SToby Isaac     PetscInt jloc = iwork[j];
3858cd1e013SToby Isaac     PetscInt diff = jloc - i;
3868cd1e013SToby Isaac 
3878cd1e013SToby Isaac     idx = idx * (n - i) + diff;
3888cd1e013SToby Isaac     /* swap (i, jloc) */
3898cd1e013SToby Isaac     work[i]     = j;
3908cd1e013SToby Isaac     work[jloc]  = icur;
3918cd1e013SToby Isaac     iwork[j]    = i;
3928cd1e013SToby Isaac     iwork[icur] = jloc;
3938cd1e013SToby Isaac     odd ^= (!!diff);
3948cd1e013SToby Isaac   }
3958cd1e013SToby Isaac   *k = idx;
3968cd1e013SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
3973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3988cd1e013SToby Isaac }
3998cd1e013SToby Isaac 
4008cd1e013SToby Isaac /*MC
401fad4db65SToby Isaac    PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k).
402fad4db65SToby Isaac    The encoding is in lexicographic order.
403fad4db65SToby Isaac 
4044165533cSJose E. Roman    Input Parameters:
405fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
406fad4db65SToby Isaac .  k - an integer in [0, n]
407fad4db65SToby Isaac -  j - an index in [0, n choose k)
408fad4db65SToby Isaac 
4094165533cSJose E. Roman    Output Parameter:
410fad4db65SToby Isaac .  subset - the jth subset of size k of the integers [0, ..., n - 1]
411fad4db65SToby Isaac 
41287497f52SBarry Smith    Note:
41387497f52SBarry Smith    Limited by arguments such that n choose k can be represented by `PetscInt`
414fad4db65SToby Isaac 
415fad4db65SToby Isaac    Level: beginner
416fad4db65SToby Isaac 
417db781477SPatrick Sanan .seealso: `PetscDTSubsetIndex()`
418fad4db65SToby Isaac M*/
419d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
420d71ae5a4SJacob Faibussowitsch {
4215f80ce2aSJacob Faibussowitsch   PetscInt Nk;
4221a989b97SToby Isaac 
4231a989b97SToby Isaac   PetscFunctionBeginHot;
4249566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
4255f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
4261a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
4271a989b97SToby Isaac     PetscInt Nminusk      = Nk - Nminuskminus;
4281a989b97SToby Isaac 
4291a989b97SToby Isaac     if (j < Nminuskminus) {
4301a989b97SToby Isaac       subset[l++] = i;
4311a989b97SToby Isaac       Nk          = Nminuskminus;
4321a989b97SToby Isaac     } else {
4331a989b97SToby Isaac       j -= Nminuskminus;
4341a989b97SToby Isaac       Nk = Nminusk;
4351a989b97SToby Isaac     }
4361a989b97SToby Isaac   }
4373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4381a989b97SToby Isaac }
4391a989b97SToby Isaac 
440fad4db65SToby Isaac /*MC
44187497f52SBarry 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.
44287497f52SBarry Smith    This is the inverse of `PetscDTEnumSubset`.
443fad4db65SToby Isaac 
4444165533cSJose E. Roman    Input Parameters:
445fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
446fad4db65SToby Isaac .  k - an integer in [0, n]
447fad4db65SToby Isaac -  subset - an ordered subset of the integers [0, ..., n - 1]
448fad4db65SToby Isaac 
4494165533cSJose E. Roman    Output Parameter:
450fad4db65SToby Isaac .  index - the rank of the subset in lexicographic order
451fad4db65SToby Isaac 
45287497f52SBarry Smith    Note:
45387497f52SBarry Smith    Limited by arguments such that n choose k can be represented by `PetscInt`
454fad4db65SToby Isaac 
455fad4db65SToby Isaac    Level: beginner
456fad4db65SToby Isaac 
457db781477SPatrick Sanan .seealso: `PetscDTEnumSubset()`
458fad4db65SToby Isaac M*/
459d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
460d71ae5a4SJacob Faibussowitsch {
4615f80ce2aSJacob Faibussowitsch   PetscInt j = 0, Nk;
4621a989b97SToby Isaac 
46328222859SToby Isaac   PetscFunctionBegin;
46428222859SToby Isaac   *index = -1;
4659566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
4665f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
4671a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
4681a989b97SToby Isaac     PetscInt Nminusk      = Nk - Nminuskminus;
4691a989b97SToby Isaac 
4701a989b97SToby Isaac     if (subset[l] == i) {
4711a989b97SToby Isaac       l++;
4721a989b97SToby Isaac       Nk = Nminuskminus;
4731a989b97SToby Isaac     } else {
4741a989b97SToby Isaac       j += Nminuskminus;
4751a989b97SToby Isaac       Nk = Nminusk;
4761a989b97SToby Isaac     }
4771a989b97SToby Isaac   }
4781a989b97SToby Isaac   *index = j;
4793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4801a989b97SToby Isaac }
4811a989b97SToby Isaac 
482fad4db65SToby Isaac /*MC
48328222859SToby 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.
484fad4db65SToby Isaac 
4854165533cSJose E. Roman    Input Parameters:
486fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
487fad4db65SToby Isaac .  k - an integer in [0, n]
488fad4db65SToby Isaac -  j - an index in [0, n choose k)
489fad4db65SToby Isaac 
4904165533cSJose E. Roman    Output Parameters:
491fad4db65SToby Isaac +  perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set.
49228222859SToby Isaac -  isOdd - if not NULL, return whether perm is an even or odd permutation.
493fad4db65SToby Isaac 
49487497f52SBarry Smith    Note:
49587497f52SBarry Smith    Limited by arguments such that n choose k can be represented by `PetscInt`
496fad4db65SToby Isaac 
497fad4db65SToby Isaac    Level: beginner
498fad4db65SToby Isaac 
499db781477SPatrick Sanan .seealso: `PetscDTEnumSubset()`, `PetscDTSubsetIndex()`
500fad4db65SToby Isaac M*/
501d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd)
502d71ae5a4SJacob Faibussowitsch {
5035f80ce2aSJacob Faibussowitsch   PetscInt  i, l, m, Nk, odd = 0;
5045f80ce2aSJacob Faibussowitsch   PetscInt *subcomp = perm + k;
5051a989b97SToby Isaac 
50628222859SToby Isaac   PetscFunctionBegin;
50728222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
5089566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
5091a989b97SToby Isaac   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
5101a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
5111a989b97SToby Isaac     PetscInt Nminusk      = Nk - Nminuskminus;
5121a989b97SToby Isaac 
5131a989b97SToby Isaac     if (j < Nminuskminus) {
514fad4db65SToby Isaac       perm[l++] = i;
5151a989b97SToby Isaac       Nk        = Nminuskminus;
5161a989b97SToby Isaac     } else {
5171a989b97SToby Isaac       subcomp[m++] = i;
5181a989b97SToby Isaac       j -= Nminuskminus;
5191a989b97SToby Isaac       odd ^= ((k - l) & 1);
5201a989b97SToby Isaac       Nk = Nminusk;
5211a989b97SToby Isaac     }
5221a989b97SToby Isaac   }
5235f80ce2aSJacob Faibussowitsch   for (; i < n; i++) subcomp[m++] = i;
5241a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
5253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5261a989b97SToby Isaac }
5271a989b97SToby Isaac 
528ef0bb6c7SMatthew G. Knepley struct _p_PetscTabulation {
529a5b23f4aSJose E. Roman   PetscInt    K;    /* Indicates a k-jet, namely tabulated derivatives up to order k */
53019815104SMartin Diehl   PetscInt    Nr;   /* The number of tabulation replicas (often 1) */
531ef0bb6c7SMatthew G. Knepley   PetscInt    Np;   /* The number of tabulation points in a replica */
532ef0bb6c7SMatthew G. Knepley   PetscInt    Nb;   /* The number of functions tabulated */
533ef0bb6c7SMatthew G. Knepley   PetscInt    Nc;   /* The number of function components */
534ef0bb6c7SMatthew G. Knepley   PetscInt    cdim; /* The coordinate dimension */
535ef0bb6c7SMatthew G. Knepley   PetscReal **T;    /* The tabulation T[K] of functions and their derivatives
536ef0bb6c7SMatthew G. Knepley                        T[0] = B[Nr*Np][Nb][Nc]:             The basis function values at quadrature points
537ef0bb6c7SMatthew G. Knepley                        T[1] = D[Nr*Np][Nb][Nc][cdim]:       The basis function derivatives at quadrature points
538ef0bb6c7SMatthew G. Knepley                        T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */
539ef0bb6c7SMatthew G. Knepley };
540ef0bb6c7SMatthew G. Knepley typedef struct _p_PetscTabulation *PetscTabulation;
541ef0bb6c7SMatthew G. Knepley 
542d6685f55SMatthew G. Knepley typedef PetscErrorCode (*PetscProbFunc)(const PetscReal[], const PetscReal[], PetscReal[]);
543d6685f55SMatthew G. Knepley 
5449371c9d4SSatish Balay typedef enum {
5459371c9d4SSatish Balay   DTPROB_DENSITY_CONSTANT,
5469371c9d4SSatish Balay   DTPROB_DENSITY_GAUSSIAN,
5479371c9d4SSatish Balay   DTPROB_DENSITY_MAXWELL_BOLTZMANN,
5489371c9d4SSatish Balay   DTPROB_NUM_DENSITY
5499371c9d4SSatish Balay } DTProbDensityType;
550d6685f55SMatthew G. Knepley PETSC_EXTERN const char *const DTProbDensityTypes[];
551d6685f55SMatthew G. Knepley 
552d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
553d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
554d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
555d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
556d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
557d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
558d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
559d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
560d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
561d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
562d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
563ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]);
564ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]);
565d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
566d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
567d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
568ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
569ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
570ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
571ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
572ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
573ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
574d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFunc *, PetscProbFunc *, PetscProbFunc *);
575d6685f55SMatthew G. Knepley 
576d6685f55SMatthew G. Knepley #include <petscvec.h>
577d6685f55SMatthew G. Knepley 
578d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFunc, PetscReal *);
579d6685f55SMatthew G. Knepley 
58037045ce4SJed Brown #endif
581