xref: /petsc/include/petscdt.h (revision d3c69ad0bb09681fda0ed612cc8275184ea14744)
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 
9ac09b921SBarry Smith /* SUBMANSEC = DT */
10ac09b921SBarry Smith 
112cd22861SMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCQUADRATURE_CLASSID;
122cd22861SMatthew G. Knepley 
1321454ff5SMatthew G. Knepley /*S
1421454ff5SMatthew G. Knepley   PetscQuadrature - Quadrature rule for integration.
1521454ff5SMatthew G. Knepley 
16329bbf4eSMatthew G. Knepley   Level: beginner
1721454ff5SMatthew G. Knepley 
18db781477SPatrick Sanan .seealso: `PetscQuadratureCreate()`, `PetscQuadratureDestroy()`
1921454ff5SMatthew G. Knepley S*/
2021454ff5SMatthew G. Knepley typedef struct _p_PetscQuadrature *PetscQuadrature;
2121454ff5SMatthew G. Knepley 
228272889dSSatish Balay /*E
23916e780bShannah_mairs   PetscGaussLobattoLegendreCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights
248272889dSSatish Balay 
258272889dSSatish Balay   Level: intermediate
268272889dSSatish Balay 
27f2e8fe4dShannah_mairs $  PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA - compute the nodes via linear algebra
28d410ae54Shannah_mairs $  PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON - compute the nodes by solving a nonlinear equation with Newton's method
298272889dSSatish Balay 
308272889dSSatish Balay E*/
31f2e8fe4dShannah_mairs typedef enum {PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON} PetscGaussLobattoLegendreCreateType;
328272889dSSatish Balay 
33d4afb720SToby Isaac /*E
34d4afb720SToby Isaac   PetscDTNodeType - A description of strategies for generating nodes (both
35d4afb720SToby Isaac   quadrature nodes and nodes for Lagrange polynomials)
36d4afb720SToby Isaac 
37d4afb720SToby Isaac   Level: intermediate
38d4afb720SToby Isaac 
39d4afb720SToby Isaac $  PETSCDTNODES_DEFAULT - Nodes chosen by PETSc
40d4afb720SToby Isaac $  PETSCDTNODES_GAUSSJACOBI - Nodes at either Gauss-Jacobi or Gauss-Lobatto-Jacobi quadrature points
41d4afb720SToby Isaac $  PETSCDTNODES_EQUISPACED - Nodes equispaced either including the endpoints or excluding them
42d4afb720SToby Isaac $  PETSCDTNODES_TANHSINH - Nodes at Tanh-Sinh quadrature points
43d4afb720SToby Isaac 
44d4afb720SToby Isaac   Note: a PetscDTNodeType can be paired with a PetscBool to indicate whether
45d4afb720SToby Isaac   the nodes include endpoints or not, and in the case of PETSCDT_GAUSSJACOBI
46d4afb720SToby Isaac   with exponents for the weight function.
47d4afb720SToby Isaac 
48d4afb720SToby Isaac E*/
49d4afb720SToby Isaac typedef enum {PETSCDTNODES_DEFAULT=-1, PETSCDTNODES_GAUSSJACOBI, PETSCDTNODES_EQUISPACED, PETSCDTNODES_TANHSINH} PetscDTNodeType;
50d4afb720SToby Isaac 
51*d3c69ad0SToby Isaac PETSC_EXTERN const char *const*const PetscDTNodeTypes;
52*d3c69ad0SToby Isaac 
53*d3c69ad0SToby Isaac /*E
54*d3c69ad0SToby Isaac   PetscDTSimplexQuadratureType - A description of classes of quadrature rules for simplices
55*d3c69ad0SToby Isaac 
56*d3c69ad0SToby Isaac   Level: intermediate
57*d3c69ad0SToby Isaac 
58*d3c69ad0SToby Isaac $  PETSCDTSIMPLEXQUAD_DEFAULT - Quadrature rule chosen by PETSc
59*d3c69ad0SToby Isaac $  PETSCDTSIMPLEXQUAD_CONIC   - Quadrature rules constructed as
60*d3c69ad0SToby Isaac                                 conically-warped tensor products of 1D
61*d3c69ad0SToby Isaac                                 Gauss-Jacobi quadrature rules.  These are
62*d3c69ad0SToby Isaac                                 explicitly computable in any dimension for any
63*d3c69ad0SToby Isaac                                 degree, and the tensor-product structure can be
64*d3c69ad0SToby Isaac                                 exploited by sum-factorization methods, but
65*d3c69ad0SToby Isaac                                 they are not efficient in terms of nodes per
66*d3c69ad0SToby Isaac                                 polynomial degree.
67*d3c69ad0SToby Isaac $  PETSCDTSIMPLEXQUAD_MINSYM  - Quadrature rules that are fully symmetric
68*d3c69ad0SToby Isaac                                 (symmetries of the simplex preserve the nodes
69*d3c69ad0SToby Isaac                                 and weights) with minimal (or near minimal)
70*d3c69ad0SToby Isaac                                 number of nodes.  In dimensions higher than 1
71*d3c69ad0SToby Isaac                                 these are not simple to compute, so lookup
72*d3c69ad0SToby Isaac                                 tables are used.
73*d3c69ad0SToby Isaac 
74*d3c69ad0SToby Isaac .seealso: `PetscDTSimplexQuadrature()`
75*d3c69ad0SToby Isaac E*/
76*d3c69ad0SToby Isaac typedef enum {PETSCDTSIMPLEXQUAD_DEFAULT=-1, PETSCDTSIMPLEXQUAD_CONIC=0, PETSCDTSIMPLEXQUAD_MINSYM} PetscDTSimplexQuadratureType;
77*d3c69ad0SToby Isaac 
78*d3c69ad0SToby Isaac PETSC_EXTERN const char *const*const PetscDTSimplexQuadratureTypes;
79d4afb720SToby Isaac 
8021454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *);
81c9638911SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *);
82bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt*);
83bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt);
84a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt*);
85a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt);
864f9ab2b4SJed Brown PETSC_EXTERN PetscErrorCode PetscQuadratureEqual(PetscQuadrature, PetscQuadrature, PetscBool*);
87a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt*, PetscInt*, PetscInt*, const PetscReal *[], const PetscReal *[]);
88a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal [], const PetscReal []);
8921454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer);
9021454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *);
91a0845e3aSMatthew G. Knepley 
922df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTensorQuadratureCreate(PetscQuadrature, PetscQuadrature, PetscQuadrature *);
9389710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *);
9489710940SMatthew G. Knepley 
95907761f8SToby Isaac PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *);
96907761f8SToby Isaac 
9737045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*);
98fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiNorm(PetscReal,PetscReal,PetscInt,PetscReal *);
9994e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt,PetscReal,PetscReal,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*);
100fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEvalJet(PetscReal,PetscReal,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscReal[]);
101fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPKDEvalJet(PetscInt,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscReal[]);
102d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedSize(PetscInt,PetscInt,PetscInt,PetscInt*);
103d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscInt,PetscReal[]);
10437045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt,PetscReal,PetscReal,PetscReal*,PetscReal*);
10594e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*,PetscReal*);
10694e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*,PetscReal*);
107916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt,PetscGaussLobattoLegendreCreateType,PetscReal*,PetscReal*);
108194825f6SJed Brown PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt,PetscInt,const PetscReal*,PetscInt,const PetscReal*,PetscReal*);
109a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
110e6a796c3SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
111*d3c69ad0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTSimplexQuadrature(PetscInt,PetscInt,PetscDTSimplexQuadratureType,PetscQuadrature*);
11237045ce4SJed Brown 
113b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
114d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
115d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
116b3c0f97bSTom Klotz 
117916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *);
118916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
119916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
120916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
121916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
122916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
123916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
124916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
125916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
126916e780bShannah_mairs 
1271a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1281a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1291a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
1301a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
1311a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
1321a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1331a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
134dda711d0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
1351a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
1361a989b97SToby Isaac 
137d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt,PetscInt,const PetscInt[],PetscInt*);
138d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt,PetscInt,PetscInt,PetscInt[]);
139fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGradedOrderToIndex(PetscInt,const PetscInt[],PetscInt*);
140fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToGradedOrder(PetscInt,PetscInt,PetscInt[]);
141d4afb720SToby Isaac 
142fad4db65SToby Isaac #if defined(PETSC_USE_64BIT_INDICES)
143fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 20
144fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX  61
145fad4db65SToby Isaac #else
146fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 12
147fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX  29
148fad4db65SToby Isaac #endif
149fad4db65SToby Isaac 
150fad4db65SToby Isaac /*MC
151fad4db65SToby Isaac    PetscDTFactorial - Approximate n! as a real number
152fad4db65SToby Isaac 
1534165533cSJose E. Roman    Input Parameter:
154fad4db65SToby Isaac .  n - a non-negative integer
155fad4db65SToby Isaac 
1564165533cSJose E. Roman    Output Parameter:
157fad4db65SToby Isaac .  factorial - n!
158fad4db65SToby Isaac 
159fad4db65SToby Isaac    Level: beginner
160fad4db65SToby Isaac M*/
1619fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial)
162fad4db65SToby Isaac {
163fad4db65SToby Isaac   PetscReal f = 1.0;
164fad4db65SToby Isaac 
165fad4db65SToby Isaac   PetscFunctionBegin;
166e2ab39ccSLisandro Dalcin   *factorial = -1.0;
16763a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %" PetscInt_FMT, n);
1685f80ce2aSJacob Faibussowitsch   for (PetscInt i = 1; i < n+1; ++i) f *= (PetscReal)i;
169fad4db65SToby Isaac   *factorial = f;
170fad4db65SToby Isaac   PetscFunctionReturn(0);
171fad4db65SToby Isaac }
172fad4db65SToby Isaac 
173fad4db65SToby Isaac /*MC
174fad4db65SToby Isaac    PetscDTFactorialInt - Compute n! as an integer
175fad4db65SToby Isaac 
1764165533cSJose E. Roman    Input Parameter:
177fad4db65SToby Isaac .  n - a non-negative integer
178fad4db65SToby Isaac 
1794165533cSJose E. Roman    Output Parameter:
180fad4db65SToby Isaac .  factorial - n!
181fad4db65SToby Isaac 
182fad4db65SToby Isaac    Level: beginner
183fad4db65SToby Isaac 
184fad4db65SToby 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.
185fad4db65SToby Isaac M*/
1869fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial)
187fad4db65SToby Isaac {
188fad4db65SToby Isaac   PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600};
189fad4db65SToby Isaac 
19028222859SToby Isaac   PetscFunctionBegin;
19128222859SToby Isaac   *factorial = -1;
19263a3b9bcSJacob 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);
193fad4db65SToby Isaac   if (n <= 12) {
194fad4db65SToby Isaac     *factorial = facLookup[n];
195fad4db65SToby Isaac   } else {
196fad4db65SToby Isaac     PetscInt f = facLookup[12];
197fad4db65SToby Isaac     PetscInt i;
198fad4db65SToby Isaac 
199fad4db65SToby Isaac     for (i = 13; i < n+1; ++i) f *= i;
200fad4db65SToby Isaac     *factorial = f;
201fad4db65SToby Isaac   }
202fad4db65SToby Isaac   PetscFunctionReturn(0);
203fad4db65SToby Isaac }
204fad4db65SToby Isaac 
205fad4db65SToby Isaac /*MC
206fad4db65SToby Isaac    PetscDTBinomial - Approximate the binomial coefficient "n choose k"
207fad4db65SToby Isaac 
2084165533cSJose E. Roman    Input Parameters:
209fad4db65SToby Isaac +  n - a non-negative integer
210fad4db65SToby Isaac -  k - an integer between 0 and n, inclusive
211fad4db65SToby Isaac 
2124165533cSJose E. Roman    Output Parameter:
213fad4db65SToby Isaac .  binomial - approximation of the binomial coefficient n choose k
214fad4db65SToby Isaac 
215fad4db65SToby Isaac    Level: beginner
216fad4db65SToby Isaac M*/
2179fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
2181a989b97SToby Isaac {
2191a989b97SToby Isaac   PetscFunctionBeginHot;
220e2ab39ccSLisandro Dalcin   *binomial = -1.0;
22163a3b9bcSJacob 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);
2221a989b97SToby Isaac   if (n <= 3) {
2231a989b97SToby Isaac     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};
2241a989b97SToby Isaac 
225e2ab39ccSLisandro Dalcin     *binomial = (PetscReal)binomLookup[n][k];
2261a989b97SToby Isaac   } else {
227e2ab39ccSLisandro Dalcin     PetscReal binom = 1.0;
2281a989b97SToby Isaac 
2291a989b97SToby Isaac     k = PetscMin(k, n - k);
2305f80ce2aSJacob Faibussowitsch     for (PetscInt i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1);
2311a989b97SToby Isaac     *binomial = binom;
2321a989b97SToby Isaac   }
2331a989b97SToby Isaac   PetscFunctionReturn(0);
2341a989b97SToby Isaac }
2351a989b97SToby Isaac 
236fad4db65SToby Isaac /*MC
237fad4db65SToby Isaac    PetscDTBinomialInt - Compute the binomial coefficient "n choose k"
238fad4db65SToby Isaac 
23997bb3fdcSJose E. Roman    Input Parameters:
240fad4db65SToby Isaac +  n - a non-negative integer
241fad4db65SToby Isaac -  k - an integer between 0 and n, inclusive
242fad4db65SToby Isaac 
24397bb3fdcSJose E. Roman    Output Parameter:
244fad4db65SToby Isaac .  binomial - the binomial coefficient n choose k
245fad4db65SToby Isaac 
246fad4db65SToby Isaac    Note: this is limited by integers that can be represented by PetscInt
247fad4db65SToby Isaac 
248fad4db65SToby Isaac    Level: beginner
249fad4db65SToby Isaac M*/
2509fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial)
251fad4db65SToby Isaac {
25228222859SToby Isaac   PetscInt bin;
25328222859SToby Isaac 
25428222859SToby Isaac   PetscFunctionBegin;
25528222859SToby Isaac   *binomial = -1;
25663a3b9bcSJacob 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);
25763a3b9bcSJacob 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);
258fad4db65SToby Isaac   if (n <= 3) {
259fad4db65SToby Isaac     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};
260fad4db65SToby Isaac 
26128222859SToby Isaac     bin = binomLookup[n][k];
262fad4db65SToby Isaac   } else {
263fad4db65SToby Isaac     PetscInt binom = 1;
264fad4db65SToby Isaac 
265fad4db65SToby Isaac     k = PetscMin(k, n - k);
2665f80ce2aSJacob Faibussowitsch     for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
26728222859SToby Isaac     bin = binom;
268fad4db65SToby Isaac   }
26928222859SToby Isaac   *binomial = bin;
270fad4db65SToby Isaac   PetscFunctionReturn(0);
271fad4db65SToby Isaac }
272fad4db65SToby Isaac 
273fad4db65SToby Isaac /*MC
274fad4db65SToby Isaac    PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps.
275fad4db65SToby Isaac 
276fad4db65SToby Isaac    A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation,
277fad4db65SToby Isaac    by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in
27828222859SToby Isaac    some position j >= i.  This swap is encoded as the difference (j - i).  The difference d_i at step i is less than
27928222859SToby Isaac    (n - i).  This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number
280fad4db65SToby Isaac    (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}.
281fad4db65SToby Isaac 
2824165533cSJose E. Roman    Input Parameters:
283fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
2848cd1e013SToby Isaac -  k - an integer in [0, n!)
285fad4db65SToby Isaac 
2864165533cSJose E. Roman    Output Parameters:
287fad4db65SToby Isaac +  perm - the permuted list of the integers [0, ..., n-1]
2888cd1e013SToby Isaac -  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.
289fad4db65SToby Isaac 
290fad4db65SToby 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.
291fad4db65SToby Isaac 
292fad4db65SToby Isaac    Level: beginner
293fad4db65SToby Isaac M*/
2949fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd)
2951a989b97SToby Isaac {
2961a989b97SToby Isaac   PetscInt  odd = 0;
2971a989b97SToby Isaac   PetscInt  i;
298fad4db65SToby Isaac   PetscInt  work[PETSC_FACTORIAL_MAX];
299fad4db65SToby Isaac   PetscInt *w;
3001a989b97SToby Isaac 
30128222859SToby Isaac   PetscFunctionBegin;
30228222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
30363a3b9bcSJacob 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);
304fad4db65SToby Isaac   w = &work[n - 2];
3051a989b97SToby Isaac   for (i = 2; i <= n; i++) {
3061a989b97SToby Isaac     *(w--) = k % i;
3071a989b97SToby Isaac     k /= i;
3081a989b97SToby Isaac   }
3091a989b97SToby Isaac   for (i = 0; i < n; i++) perm[i] = i;
3101a989b97SToby Isaac   for (i = 0; i < n - 1; i++) {
3111a989b97SToby Isaac     PetscInt s = work[i];
3121a989b97SToby Isaac     PetscInt swap = perm[i];
3131a989b97SToby Isaac 
3141a989b97SToby Isaac     perm[i] = perm[i + s];
3151a989b97SToby Isaac     perm[i + s] = swap;
3161a989b97SToby Isaac     odd ^= (!!s);
3171a989b97SToby Isaac   }
3181a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
3191a989b97SToby Isaac   PetscFunctionReturn(0);
3201a989b97SToby Isaac }
3211a989b97SToby Isaac 
322fad4db65SToby Isaac /*MC
3238cd1e013SToby Isaac    PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!).  This inverts PetscDTEnumPerm.
3248cd1e013SToby Isaac 
3254165533cSJose E. Roman    Input Parameters:
3268cd1e013SToby Isaac +  n - a non-negative integer (see note about limits below)
3278cd1e013SToby Isaac -  perm - the permuted list of the integers [0, ..., n-1]
3288cd1e013SToby Isaac 
3294165533cSJose E. Roman    Output Parameters:
3308cd1e013SToby Isaac +  k - an integer in [0, n!)
331f0fc11ceSJed Brown -  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.
3328cd1e013SToby Isaac 
3338cd1e013SToby 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.
3348cd1e013SToby Isaac 
3358cd1e013SToby Isaac    Level: beginner
3368cd1e013SToby Isaac M*/
3379fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd)
3388cd1e013SToby Isaac {
3398cd1e013SToby Isaac   PetscInt  odd = 0;
3408cd1e013SToby Isaac   PetscInt  i, idx;
3418cd1e013SToby Isaac   PetscInt  work[PETSC_FACTORIAL_MAX];
3428cd1e013SToby Isaac   PetscInt  iwork[PETSC_FACTORIAL_MAX];
3438cd1e013SToby Isaac 
3448cd1e013SToby Isaac   PetscFunctionBeginHot;
34528222859SToby Isaac   *k = -1;
34628222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
34763a3b9bcSJacob 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);
3488cd1e013SToby Isaac   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
3498cd1e013SToby Isaac   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
3508cd1e013SToby Isaac   for (idx = 0, i = 0; i < n - 1; i++) {
3518cd1e013SToby Isaac     PetscInt j = perm[i];
3528cd1e013SToby Isaac     PetscInt icur = work[i];
3538cd1e013SToby Isaac     PetscInt jloc = iwork[j];
3548cd1e013SToby Isaac     PetscInt diff = jloc - i;
3558cd1e013SToby Isaac 
3568cd1e013SToby Isaac     idx = idx * (n - i) + diff;
3578cd1e013SToby Isaac     /* swap (i, jloc) */
3588cd1e013SToby Isaac     work[i] = j;
3598cd1e013SToby Isaac     work[jloc] = icur;
3608cd1e013SToby Isaac     iwork[j] = i;
3618cd1e013SToby Isaac     iwork[icur] = jloc;
3628cd1e013SToby Isaac     odd ^= (!!diff);
3638cd1e013SToby Isaac   }
3648cd1e013SToby Isaac   *k = idx;
3658cd1e013SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
3668cd1e013SToby Isaac   PetscFunctionReturn(0);
3678cd1e013SToby Isaac }
3688cd1e013SToby Isaac 
3698cd1e013SToby Isaac /*MC
370fad4db65SToby Isaac    PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k).
371fad4db65SToby Isaac    The encoding is in lexicographic order.
372fad4db65SToby Isaac 
3734165533cSJose E. Roman    Input Parameters:
374fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
375fad4db65SToby Isaac .  k - an integer in [0, n]
376fad4db65SToby Isaac -  j - an index in [0, n choose k)
377fad4db65SToby Isaac 
3784165533cSJose E. Roman    Output Parameter:
379fad4db65SToby Isaac .  subset - the jth subset of size k of the integers [0, ..., n - 1]
380fad4db65SToby Isaac 
381fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
382fad4db65SToby Isaac 
383fad4db65SToby Isaac    Level: beginner
384fad4db65SToby Isaac 
385db781477SPatrick Sanan .seealso: `PetscDTSubsetIndex()`
386fad4db65SToby Isaac M*/
3879fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
3881a989b97SToby Isaac {
3895f80ce2aSJacob Faibussowitsch   PetscInt Nk;
3901a989b97SToby Isaac 
3911a989b97SToby Isaac   PetscFunctionBeginHot;
3929566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
3935f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
3941a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
3951a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
3961a989b97SToby Isaac 
3971a989b97SToby Isaac     if (j < Nminuskminus) {
3981a989b97SToby Isaac       subset[l++] = i;
3991a989b97SToby Isaac       Nk = Nminuskminus;
4001a989b97SToby Isaac     } else {
4011a989b97SToby Isaac       j -= Nminuskminus;
4021a989b97SToby Isaac       Nk = Nminusk;
4031a989b97SToby Isaac     }
4041a989b97SToby Isaac   }
4051a989b97SToby Isaac   PetscFunctionReturn(0);
4061a989b97SToby Isaac }
4071a989b97SToby Isaac 
408fad4db65SToby Isaac /*MC
409fad4db65SToby 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.
410fad4db65SToby Isaac 
4114165533cSJose E. Roman    Input Parameters:
412fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
413fad4db65SToby Isaac .  k - an integer in [0, n]
414fad4db65SToby Isaac -  subset - an ordered subset of the integers [0, ..., n - 1]
415fad4db65SToby Isaac 
4164165533cSJose E. Roman    Output Parameter:
417fad4db65SToby Isaac .  index - the rank of the subset in lexicographic order
418fad4db65SToby Isaac 
419fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
420fad4db65SToby Isaac 
421fad4db65SToby Isaac    Level: beginner
422fad4db65SToby Isaac 
423db781477SPatrick Sanan .seealso: `PetscDTEnumSubset()`
424fad4db65SToby Isaac M*/
4259fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
4261a989b97SToby Isaac {
4275f80ce2aSJacob Faibussowitsch   PetscInt j = 0, Nk;
4281a989b97SToby Isaac 
42928222859SToby Isaac   PetscFunctionBegin;
43028222859SToby Isaac   *index = -1;
4319566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
4325f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
4331a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
4341a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
4351a989b97SToby Isaac 
4361a989b97SToby Isaac     if (subset[l] == i) {
4371a989b97SToby Isaac       l++;
4381a989b97SToby Isaac       Nk = Nminuskminus;
4391a989b97SToby Isaac     } else {
4401a989b97SToby Isaac       j += Nminuskminus;
4411a989b97SToby Isaac       Nk = Nminusk;
4421a989b97SToby Isaac     }
4431a989b97SToby Isaac   }
4441a989b97SToby Isaac   *index = j;
4451a989b97SToby Isaac   PetscFunctionReturn(0);
4461a989b97SToby Isaac }
4471a989b97SToby Isaac 
448fad4db65SToby Isaac /*MC
44928222859SToby 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.
450fad4db65SToby Isaac 
4514165533cSJose E. Roman    Input Parameters:
452fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
453fad4db65SToby Isaac .  k - an integer in [0, n]
454fad4db65SToby Isaac -  j - an index in [0, n choose k)
455fad4db65SToby Isaac 
4564165533cSJose E. Roman    Output Parameters:
457fad4db65SToby Isaac +  perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set.
45828222859SToby Isaac -  isOdd - if not NULL, return whether perm is an even or odd permutation.
459fad4db65SToby Isaac 
460fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
461fad4db65SToby Isaac 
462fad4db65SToby Isaac    Level: beginner
463fad4db65SToby Isaac 
464db781477SPatrick Sanan .seealso: `PetscDTEnumSubset()`, `PetscDTSubsetIndex()`
465fad4db65SToby Isaac M*/
4669fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd)
4671a989b97SToby Isaac {
4685f80ce2aSJacob Faibussowitsch   PetscInt i, l, m, Nk, odd = 0;
4695f80ce2aSJacob Faibussowitsch   PetscInt *subcomp = perm+k;
4701a989b97SToby Isaac 
47128222859SToby Isaac   PetscFunctionBegin;
47228222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
4739566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
4741a989b97SToby Isaac   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
4751a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
4761a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
4771a989b97SToby Isaac 
4781a989b97SToby Isaac     if (j < Nminuskminus) {
479fad4db65SToby Isaac       perm[l++] = i;
4801a989b97SToby Isaac       Nk = Nminuskminus;
4811a989b97SToby Isaac     } else {
4821a989b97SToby Isaac       subcomp[m++] = i;
4831a989b97SToby Isaac       j -= Nminuskminus;
4841a989b97SToby Isaac       odd ^= ((k - l) & 1);
4851a989b97SToby Isaac       Nk = Nminusk;
4861a989b97SToby Isaac     }
4871a989b97SToby Isaac   }
4885f80ce2aSJacob Faibussowitsch   for (; i < n; i++) subcomp[m++] = i;
4891a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
4901a989b97SToby Isaac   PetscFunctionReturn(0);
4911a989b97SToby Isaac }
4921a989b97SToby Isaac 
493ef0bb6c7SMatthew G. Knepley struct _p_PetscTabulation {
494a5b23f4aSJose E. Roman   PetscInt    K;    /* Indicates a k-jet, namely tabulated derivatives up to order k */
49519815104SMartin Diehl   PetscInt    Nr;   /* The number of tabulation replicas (often 1) */
496ef0bb6c7SMatthew G. Knepley   PetscInt    Np;   /* The number of tabulation points in a replica */
497ef0bb6c7SMatthew G. Knepley   PetscInt    Nb;   /* The number of functions tabulated */
498ef0bb6c7SMatthew G. Knepley   PetscInt    Nc;   /* The number of function components */
499ef0bb6c7SMatthew G. Knepley   PetscInt    cdim; /* The coordinate dimension */
500ef0bb6c7SMatthew G. Knepley   PetscReal **T;    /* The tabulation T[K] of functions and their derivatives
501ef0bb6c7SMatthew G. Knepley                        T[0] = B[Nr*Np][Nb][Nc]:             The basis function values at quadrature points
502ef0bb6c7SMatthew G. Knepley                        T[1] = D[Nr*Np][Nb][Nc][cdim]:       The basis function derivatives at quadrature points
503ef0bb6c7SMatthew G. Knepley                        T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */
504ef0bb6c7SMatthew G. Knepley };
505ef0bb6c7SMatthew G. Knepley typedef struct _p_PetscTabulation *PetscTabulation;
506ef0bb6c7SMatthew G. Knepley 
507d6685f55SMatthew G. Knepley typedef PetscErrorCode (*PetscProbFunc)(const PetscReal[], const PetscReal[], PetscReal[]);
508d6685f55SMatthew G. Knepley 
509d6685f55SMatthew G. Knepley typedef enum {DTPROB_DENSITY_CONSTANT, DTPROB_DENSITY_GAUSSIAN, DTPROB_DENSITY_MAXWELL_BOLTZMANN, DTPROB_NUM_DENSITY} DTProbDensityType;
510d6685f55SMatthew G. Knepley PETSC_EXTERN const char * const DTProbDensityTypes[];
511d6685f55SMatthew G. Knepley 
512d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
513d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
514d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
515d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
516d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
517d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
518d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
519d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
520d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
521d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
522d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
523ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]);
524ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]);
525d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
526d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
527d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
528ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
529ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
530ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
531ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
532ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
533ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
534d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFunc *, PetscProbFunc *, PetscProbFunc *);
535d6685f55SMatthew G. Knepley 
536d6685f55SMatthew G. Knepley #include <petscvec.h>
537d6685f55SMatthew G. Knepley 
538d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFunc, PetscReal *);
539d6685f55SMatthew G. Knepley 
54037045ce4SJed Brown #endif
541