137045ce4SJed Brown /* 237045ce4SJed Brown Common tools for constructing discretizations 337045ce4SJed Brown */ 4*6524c165SJacob Faibussowitsch #ifndef 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 2787497f52SBarry Smith $ `PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA` - compute the nodes via linear algebra 2887497f52SBarry Smith $ `PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON` - compute the nodes by solving a nonlinear equation with Newton's method 298272889dSSatish Balay 308272889dSSatish Balay E*/ 319371c9d4SSatish Balay typedef enum { 329371c9d4SSatish Balay PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA, 339371c9d4SSatish Balay PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON 349371c9d4SSatish Balay } PetscGaussLobattoLegendreCreateType; 358272889dSSatish Balay 36d4afb720SToby Isaac /*E 37d4afb720SToby Isaac PetscDTNodeType - A description of strategies for generating nodes (both 38d4afb720SToby Isaac quadrature nodes and nodes for Lagrange polynomials) 39d4afb720SToby Isaac 40d4afb720SToby Isaac Level: intermediate 41d4afb720SToby Isaac 4287497f52SBarry Smith $ `PETSCDTNODES_DEFAULT` - Nodes chosen by PETSc 4387497f52SBarry Smith $ `PETSCDTNODES_GAUSSJACOBI` - Nodes at either Gauss-Jacobi or Gauss-Lobatto-Jacobi quadrature points 4487497f52SBarry Smith $ `PETSCDTNODES_EQUISPACED` - Nodes equispaced either including the endpoints or excluding them 4587497f52SBarry Smith $ `PETSCDTNODES_TANHSINH` - Nodes at Tanh-Sinh quadrature points 46d4afb720SToby Isaac 4787497f52SBarry Smith Note: 4887497f52SBarry Smith A `PetscDTNodeType` can be paired with a `PetscBool` to indicate whether 4987497f52SBarry Smith the nodes include endpoints or not, and in the case of `PETSCDT_GAUSSJACOBI` 50d4afb720SToby Isaac with exponents for the weight function. 51d4afb720SToby Isaac 52d4afb720SToby Isaac E*/ 539371c9d4SSatish Balay typedef enum { 549371c9d4SSatish Balay PETSCDTNODES_DEFAULT = -1, 559371c9d4SSatish Balay PETSCDTNODES_GAUSSJACOBI, 569371c9d4SSatish Balay PETSCDTNODES_EQUISPACED, 579371c9d4SSatish Balay PETSCDTNODES_TANHSINH 589371c9d4SSatish Balay } PetscDTNodeType; 59d4afb720SToby Isaac 60d3c69ad0SToby Isaac PETSC_EXTERN const char *const *const PetscDTNodeTypes; 61d3c69ad0SToby Isaac 62d3c69ad0SToby Isaac /*E 63d3c69ad0SToby Isaac PetscDTSimplexQuadratureType - A description of classes of quadrature rules for simplices 64d3c69ad0SToby Isaac 65d3c69ad0SToby Isaac Level: intermediate 66d3c69ad0SToby Isaac 6787497f52SBarry Smith $ `PETSCDTSIMPLEXQUAD_DEFAULT` - Quadrature rule chosen by PETSc 6887497f52SBarry Smith $ `PETSCDTSIMPLEXQUAD_CONIC` - Quadrature rules constructed as 69d3c69ad0SToby Isaac conically-warped tensor products of 1D 70d3c69ad0SToby Isaac Gauss-Jacobi quadrature rules. These are 71d3c69ad0SToby Isaac explicitly computable in any dimension for any 72d3c69ad0SToby Isaac degree, and the tensor-product structure can be 73d3c69ad0SToby Isaac exploited by sum-factorization methods, but 74d3c69ad0SToby Isaac they are not efficient in terms of nodes per 75d3c69ad0SToby Isaac polynomial degree. 7687497f52SBarry Smith $ `PETSCDTSIMPLEXQUAD_MINSYM` - Quadrature rules that are fully symmetric 77d3c69ad0SToby Isaac (symmetries of the simplex preserve the nodes 78d3c69ad0SToby Isaac and weights) with minimal (or near minimal) 79d3c69ad0SToby Isaac number of nodes. In dimensions higher than 1 80d3c69ad0SToby Isaac these are not simple to compute, so lookup 81d3c69ad0SToby Isaac tables are used. 82d3c69ad0SToby Isaac 83d3c69ad0SToby Isaac .seealso: `PetscDTSimplexQuadrature()` 84d3c69ad0SToby Isaac E*/ 859371c9d4SSatish Balay typedef enum { 869371c9d4SSatish Balay PETSCDTSIMPLEXQUAD_DEFAULT = -1, 879371c9d4SSatish Balay PETSCDTSIMPLEXQUAD_CONIC = 0, 889371c9d4SSatish Balay PETSCDTSIMPLEXQUAD_MINSYM 899371c9d4SSatish Balay } PetscDTSimplexQuadratureType; 90d3c69ad0SToby Isaac 91d3c69ad0SToby Isaac PETSC_EXTERN const char *const *const PetscDTSimplexQuadratureTypes; 92d4afb720SToby Isaac 9321454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *); 94c9638911SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *); 95bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt *); 96bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt); 97a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt *); 98a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt); 994f9ab2b4SJed Brown PETSC_EXTERN PetscErrorCode PetscQuadratureEqual(PetscQuadrature, PetscQuadrature, PetscBool *); 100a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt *, PetscInt *, PetscInt *, const PetscReal *[], const PetscReal *[]); 101a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal[], const PetscReal[]); 10221454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer); 10321454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *); 104a0845e3aSMatthew G. Knepley 1052df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTensorQuadratureCreate(PetscQuadrature, PetscQuadrature, PetscQuadrature *); 10689710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *); 10789710940SMatthew G. Knepley 108907761f8SToby Isaac PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *); 109907761f8SToby Isaac 11037045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *); 111fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiNorm(PetscReal, PetscReal, PetscInt, PetscReal *); 11294e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt, PetscReal, PetscReal, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *); 113fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEvalJet(PetscReal, PetscReal, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]); 114fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPKDEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]); 115d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedSize(PetscInt, PetscInt, PetscInt, PetscInt *); 116d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscInt, PetscReal[]); 11737045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt, PetscReal, PetscReal, PetscReal *, PetscReal *); 11894e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *); 11994e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *); 120916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt, PetscGaussLobattoLegendreCreateType, PetscReal *, PetscReal *); 121194825f6SJed Brown PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *); 122a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *); 123e6a796c3SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *); 124d3c69ad0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTSimplexQuadrature(PetscInt, PetscInt, PetscDTSimplexQuadratureType, PetscQuadrature *); 12537045ce4SJed Brown 126b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *); 127d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *); 128d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *); 129b3c0f97bSTom Klotz 130916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *); 131916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 132916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 133916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 134916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 135916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 136916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 137916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 138916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 139916e780bShannah_mairs 1401a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 1411a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 1421a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 1431a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *); 1441a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *); 1451a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 1461a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *); 147dda711d0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]); 1481a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 1491a989b97SToby Isaac 150d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt, PetscInt, const PetscInt[], PetscInt *); 151d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt, PetscInt, PetscInt, PetscInt[]); 152fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGradedOrderToIndex(PetscInt, const PetscInt[], PetscInt *); 153fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToGradedOrder(PetscInt, PetscInt, PetscInt[]); 154d4afb720SToby Isaac 155fad4db65SToby Isaac #if defined(PETSC_USE_64BIT_INDICES) 156fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 20 157fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX 61 158fad4db65SToby Isaac #else 159fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 12 160fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX 29 161fad4db65SToby Isaac #endif 162fad4db65SToby Isaac 163fad4db65SToby Isaac /*MC 164fad4db65SToby Isaac PetscDTFactorial - Approximate n! as a real number 165fad4db65SToby Isaac 1664165533cSJose E. Roman Input Parameter: 167fad4db65SToby Isaac . n - a non-negative integer 168fad4db65SToby Isaac 1694165533cSJose E. Roman Output Parameter: 170fad4db65SToby Isaac . factorial - n! 171fad4db65SToby Isaac 172fad4db65SToby Isaac Level: beginner 173fad4db65SToby Isaac M*/ 1749371c9d4SSatish Balay static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial) { 175fad4db65SToby Isaac PetscReal f = 1.0; 176fad4db65SToby Isaac 177fad4db65SToby Isaac PetscFunctionBegin; 178e2ab39ccSLisandro Dalcin *factorial = -1.0; 17963a3b9bcSJacob Faibussowitsch PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %" PetscInt_FMT, n); 1805f80ce2aSJacob Faibussowitsch for (PetscInt i = 1; i < n + 1; ++i) f *= (PetscReal)i; 181fad4db65SToby Isaac *factorial = f; 182fad4db65SToby Isaac PetscFunctionReturn(0); 183fad4db65SToby Isaac } 184fad4db65SToby Isaac 185fad4db65SToby Isaac /*MC 186fad4db65SToby Isaac PetscDTFactorialInt - Compute n! as an integer 187fad4db65SToby Isaac 1884165533cSJose E. Roman Input Parameter: 189fad4db65SToby Isaac . n - a non-negative integer 190fad4db65SToby Isaac 1914165533cSJose E. Roman Output Parameter: 192fad4db65SToby Isaac . factorial - n! 193fad4db65SToby Isaac 194fad4db65SToby Isaac Level: beginner 195fad4db65SToby Isaac 196fad4db65SToby 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. 197fad4db65SToby Isaac M*/ 1989371c9d4SSatish Balay static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial) { 199fad4db65SToby Isaac PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600}; 200fad4db65SToby Isaac 20128222859SToby Isaac PetscFunctionBegin; 20228222859SToby Isaac *factorial = -1; 20363a3b9bcSJacob 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); 204fad4db65SToby Isaac if (n <= 12) { 205fad4db65SToby Isaac *factorial = facLookup[n]; 206fad4db65SToby Isaac } else { 207fad4db65SToby Isaac PetscInt f = facLookup[12]; 208fad4db65SToby Isaac PetscInt i; 209fad4db65SToby Isaac 210fad4db65SToby Isaac for (i = 13; i < n + 1; ++i) f *= i; 211fad4db65SToby Isaac *factorial = f; 212fad4db65SToby Isaac } 213fad4db65SToby Isaac PetscFunctionReturn(0); 214fad4db65SToby Isaac } 215fad4db65SToby Isaac 216fad4db65SToby Isaac /*MC 217fad4db65SToby Isaac PetscDTBinomial - Approximate the binomial coefficient "n choose k" 218fad4db65SToby Isaac 2194165533cSJose E. Roman Input Parameters: 220fad4db65SToby Isaac + n - a non-negative integer 221fad4db65SToby Isaac - k - an integer between 0 and n, inclusive 222fad4db65SToby Isaac 2234165533cSJose E. Roman Output Parameter: 224fad4db65SToby Isaac . binomial - approximation of the binomial coefficient n choose k 225fad4db65SToby Isaac 226fad4db65SToby Isaac Level: beginner 227fad4db65SToby Isaac M*/ 2289371c9d4SSatish Balay static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial) { 2291a989b97SToby Isaac PetscFunctionBeginHot; 230e2ab39ccSLisandro Dalcin *binomial = -1.0; 23163a3b9bcSJacob 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); 2321a989b97SToby Isaac if (n <= 3) { 2339371c9d4SSatish Balay PetscInt binomLookup[4][4] = { 2349371c9d4SSatish Balay {1, 0, 0, 0}, 2359371c9d4SSatish Balay {1, 1, 0, 0}, 2369371c9d4SSatish Balay {1, 2, 1, 0}, 2379371c9d4SSatish Balay {1, 3, 3, 1} 2389371c9d4SSatish Balay }; 2391a989b97SToby Isaac 240e2ab39ccSLisandro Dalcin *binomial = (PetscReal)binomLookup[n][k]; 2411a989b97SToby Isaac } else { 242e2ab39ccSLisandro Dalcin PetscReal binom = 1.0; 2431a989b97SToby Isaac 2441a989b97SToby Isaac k = PetscMin(k, n - k); 2455f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1); 2461a989b97SToby Isaac *binomial = binom; 2471a989b97SToby Isaac } 2481a989b97SToby Isaac PetscFunctionReturn(0); 2491a989b97SToby Isaac } 2501a989b97SToby Isaac 251fad4db65SToby Isaac /*MC 252fad4db65SToby Isaac PetscDTBinomialInt - Compute the binomial coefficient "n choose k" 253fad4db65SToby Isaac 25497bb3fdcSJose E. Roman Input Parameters: 255fad4db65SToby Isaac + n - a non-negative integer 256fad4db65SToby Isaac - k - an integer between 0 and n, inclusive 257fad4db65SToby Isaac 25897bb3fdcSJose E. Roman Output Parameter: 259fad4db65SToby Isaac . binomial - the binomial coefficient n choose k 260fad4db65SToby Isaac 26187497f52SBarry Smith Note: this is limited by integers that can be represented by `PetscInt` 262fad4db65SToby Isaac 263fad4db65SToby Isaac Level: beginner 264fad4db65SToby Isaac M*/ 2659371c9d4SSatish Balay static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial) { 26628222859SToby Isaac PetscInt bin; 26728222859SToby Isaac 26828222859SToby Isaac PetscFunctionBegin; 26928222859SToby Isaac *binomial = -1; 27063a3b9bcSJacob 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); 27163a3b9bcSJacob 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); 272fad4db65SToby Isaac if (n <= 3) { 2739371c9d4SSatish Balay PetscInt binomLookup[4][4] = { 2749371c9d4SSatish Balay {1, 0, 0, 0}, 2759371c9d4SSatish Balay {1, 1, 0, 0}, 2769371c9d4SSatish Balay {1, 2, 1, 0}, 2779371c9d4SSatish Balay {1, 3, 3, 1} 2789371c9d4SSatish Balay }; 279fad4db65SToby Isaac 28028222859SToby Isaac bin = binomLookup[n][k]; 281fad4db65SToby Isaac } else { 282fad4db65SToby Isaac PetscInt binom = 1; 283fad4db65SToby Isaac 284fad4db65SToby Isaac k = PetscMin(k, n - k); 2855f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1); 28628222859SToby Isaac bin = binom; 287fad4db65SToby Isaac } 28828222859SToby Isaac *binomial = bin; 289fad4db65SToby Isaac PetscFunctionReturn(0); 290fad4db65SToby Isaac } 291fad4db65SToby Isaac 292fad4db65SToby Isaac /*MC 293fad4db65SToby Isaac PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps. 294fad4db65SToby Isaac 295fad4db65SToby Isaac A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation, 296fad4db65SToby Isaac by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in 29728222859SToby Isaac some position j >= i. This swap is encoded as the difference (j - i). The difference d_i at step i is less than 29828222859SToby Isaac (n - i). This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number 299fad4db65SToby Isaac (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}. 300fad4db65SToby Isaac 3014165533cSJose E. Roman Input Parameters: 302fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 3038cd1e013SToby Isaac - k - an integer in [0, n!) 304fad4db65SToby Isaac 3054165533cSJose E. Roman Output Parameters: 306fad4db65SToby Isaac + perm - the permuted list of the integers [0, ..., n-1] 3078cd1e013SToby Isaac - isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps. 308fad4db65SToby Isaac 30987497f52SBarry Smith Note: 31087497f52SBarry 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. 311fad4db65SToby Isaac 312fad4db65SToby Isaac Level: beginner 313fad4db65SToby Isaac M*/ 3149371c9d4SSatish Balay static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd) { 3151a989b97SToby Isaac PetscInt odd = 0; 3161a989b97SToby Isaac PetscInt i; 317fad4db65SToby Isaac PetscInt work[PETSC_FACTORIAL_MAX]; 318fad4db65SToby Isaac PetscInt *w; 3191a989b97SToby Isaac 32028222859SToby Isaac PetscFunctionBegin; 32128222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 32263a3b9bcSJacob 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); 323fad4db65SToby Isaac w = &work[n - 2]; 3241a989b97SToby Isaac for (i = 2; i <= n; i++) { 3251a989b97SToby Isaac *(w--) = k % i; 3261a989b97SToby Isaac k /= i; 3271a989b97SToby Isaac } 3281a989b97SToby Isaac for (i = 0; i < n; i++) perm[i] = i; 3291a989b97SToby Isaac for (i = 0; i < n - 1; i++) { 3301a989b97SToby Isaac PetscInt s = work[i]; 3311a989b97SToby Isaac PetscInt swap = perm[i]; 3321a989b97SToby Isaac 3331a989b97SToby Isaac perm[i] = perm[i + s]; 3341a989b97SToby Isaac perm[i + s] = swap; 3351a989b97SToby Isaac odd ^= (!!s); 3361a989b97SToby Isaac } 3371a989b97SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 3381a989b97SToby Isaac PetscFunctionReturn(0); 3391a989b97SToby Isaac } 3401a989b97SToby Isaac 341fad4db65SToby Isaac /*MC 34287497f52SBarry Smith PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!). This inverts `PetscDTEnumPerm`. 3438cd1e013SToby Isaac 3444165533cSJose E. Roman Input Parameters: 3458cd1e013SToby Isaac + n - a non-negative integer (see note about limits below) 3468cd1e013SToby Isaac - perm - the permuted list of the integers [0, ..., n-1] 3478cd1e013SToby Isaac 3484165533cSJose E. Roman Output Parameters: 3498cd1e013SToby Isaac + k - an integer in [0, n!) 350f0fc11ceSJed Brown - isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps. 3518cd1e013SToby Isaac 35287497f52SBarry Smith Note: 35387497f52SBarry 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. 3548cd1e013SToby Isaac 3558cd1e013SToby Isaac Level: beginner 3568cd1e013SToby Isaac M*/ 3579371c9d4SSatish Balay static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd) { 3588cd1e013SToby Isaac PetscInt odd = 0; 3598cd1e013SToby Isaac PetscInt i, idx; 3608cd1e013SToby Isaac PetscInt work[PETSC_FACTORIAL_MAX]; 3618cd1e013SToby Isaac PetscInt iwork[PETSC_FACTORIAL_MAX]; 3628cd1e013SToby Isaac 3638cd1e013SToby Isaac PetscFunctionBeginHot; 36428222859SToby Isaac *k = -1; 36528222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 36663a3b9bcSJacob 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); 3678cd1e013SToby Isaac for (i = 0; i < n; i++) work[i] = i; /* partial permutation */ 3688cd1e013SToby Isaac for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */ 3698cd1e013SToby Isaac for (idx = 0, i = 0; i < n - 1; i++) { 3708cd1e013SToby Isaac PetscInt j = perm[i]; 3718cd1e013SToby Isaac PetscInt icur = work[i]; 3728cd1e013SToby Isaac PetscInt jloc = iwork[j]; 3738cd1e013SToby Isaac PetscInt diff = jloc - i; 3748cd1e013SToby Isaac 3758cd1e013SToby Isaac idx = idx * (n - i) + diff; 3768cd1e013SToby Isaac /* swap (i, jloc) */ 3778cd1e013SToby Isaac work[i] = j; 3788cd1e013SToby Isaac work[jloc] = icur; 3798cd1e013SToby Isaac iwork[j] = i; 3808cd1e013SToby Isaac iwork[icur] = jloc; 3818cd1e013SToby Isaac odd ^= (!!diff); 3828cd1e013SToby Isaac } 3838cd1e013SToby Isaac *k = idx; 3848cd1e013SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 3858cd1e013SToby Isaac PetscFunctionReturn(0); 3868cd1e013SToby Isaac } 3878cd1e013SToby Isaac 3888cd1e013SToby Isaac /*MC 389fad4db65SToby Isaac PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k). 390fad4db65SToby Isaac The encoding is in lexicographic order. 391fad4db65SToby Isaac 3924165533cSJose E. Roman Input Parameters: 393fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 394fad4db65SToby Isaac . k - an integer in [0, n] 395fad4db65SToby Isaac - j - an index in [0, n choose k) 396fad4db65SToby Isaac 3974165533cSJose E. Roman Output Parameter: 398fad4db65SToby Isaac . subset - the jth subset of size k of the integers [0, ..., n - 1] 399fad4db65SToby Isaac 40087497f52SBarry Smith Note: 40187497f52SBarry Smith Limited by arguments such that n choose k can be represented by `PetscInt` 402fad4db65SToby Isaac 403fad4db65SToby Isaac Level: beginner 404fad4db65SToby Isaac 405db781477SPatrick Sanan .seealso: `PetscDTSubsetIndex()` 406fad4db65SToby Isaac M*/ 4079371c9d4SSatish Balay static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset) { 4085f80ce2aSJacob Faibussowitsch PetscInt Nk; 4091a989b97SToby Isaac 4101a989b97SToby Isaac PetscFunctionBeginHot; 4119566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 4125f80ce2aSJacob Faibussowitsch for (PetscInt i = 0, l = 0; i < n && l < k; i++) { 4131a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 4141a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 4151a989b97SToby Isaac 4161a989b97SToby Isaac if (j < Nminuskminus) { 4171a989b97SToby Isaac subset[l++] = i; 4181a989b97SToby Isaac Nk = Nminuskminus; 4191a989b97SToby Isaac } else { 4201a989b97SToby Isaac j -= Nminuskminus; 4211a989b97SToby Isaac Nk = Nminusk; 4221a989b97SToby Isaac } 4231a989b97SToby Isaac } 4241a989b97SToby Isaac PetscFunctionReturn(0); 4251a989b97SToby Isaac } 4261a989b97SToby Isaac 427fad4db65SToby Isaac /*MC 42887497f52SBarry 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. 42987497f52SBarry Smith This is the inverse of `PetscDTEnumSubset`. 430fad4db65SToby Isaac 4314165533cSJose E. Roman Input Parameters: 432fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 433fad4db65SToby Isaac . k - an integer in [0, n] 434fad4db65SToby Isaac - subset - an ordered subset of the integers [0, ..., n - 1] 435fad4db65SToby Isaac 4364165533cSJose E. Roman Output Parameter: 437fad4db65SToby Isaac . index - the rank of the subset in lexicographic order 438fad4db65SToby Isaac 43987497f52SBarry Smith Note: 44087497f52SBarry Smith Limited by arguments such that n choose k can be represented by `PetscInt` 441fad4db65SToby Isaac 442fad4db65SToby Isaac Level: beginner 443fad4db65SToby Isaac 444db781477SPatrick Sanan .seealso: `PetscDTEnumSubset()` 445fad4db65SToby Isaac M*/ 4469371c9d4SSatish Balay static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index) { 4475f80ce2aSJacob Faibussowitsch PetscInt j = 0, Nk; 4481a989b97SToby Isaac 44928222859SToby Isaac PetscFunctionBegin; 45028222859SToby Isaac *index = -1; 4519566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 4525f80ce2aSJacob Faibussowitsch for (PetscInt i = 0, l = 0; i < n && l < k; i++) { 4531a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 4541a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 4551a989b97SToby Isaac 4561a989b97SToby Isaac if (subset[l] == i) { 4571a989b97SToby Isaac l++; 4581a989b97SToby Isaac Nk = Nminuskminus; 4591a989b97SToby Isaac } else { 4601a989b97SToby Isaac j += Nminuskminus; 4611a989b97SToby Isaac Nk = Nminusk; 4621a989b97SToby Isaac } 4631a989b97SToby Isaac } 4641a989b97SToby Isaac *index = j; 4651a989b97SToby Isaac PetscFunctionReturn(0); 4661a989b97SToby Isaac } 4671a989b97SToby Isaac 468fad4db65SToby Isaac /*MC 46928222859SToby 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. 470fad4db65SToby Isaac 4714165533cSJose E. Roman Input Parameters: 472fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 473fad4db65SToby Isaac . k - an integer in [0, n] 474fad4db65SToby Isaac - j - an index in [0, n choose k) 475fad4db65SToby Isaac 4764165533cSJose E. Roman Output Parameters: 477fad4db65SToby Isaac + perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set. 47828222859SToby Isaac - isOdd - if not NULL, return whether perm is an even or odd permutation. 479fad4db65SToby Isaac 48087497f52SBarry Smith Note: 48187497f52SBarry Smith Limited by arguments such that n choose k can be represented by `PetscInt` 482fad4db65SToby Isaac 483fad4db65SToby Isaac Level: beginner 484fad4db65SToby Isaac 485db781477SPatrick Sanan .seealso: `PetscDTEnumSubset()`, `PetscDTSubsetIndex()` 486fad4db65SToby Isaac M*/ 4879371c9d4SSatish Balay static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd) { 4885f80ce2aSJacob Faibussowitsch PetscInt i, l, m, Nk, odd = 0; 4895f80ce2aSJacob Faibussowitsch PetscInt *subcomp = perm + k; 4901a989b97SToby Isaac 49128222859SToby Isaac PetscFunctionBegin; 49228222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 4939566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 4941a989b97SToby Isaac for (i = 0, l = 0, m = 0; i < n && l < k; i++) { 4951a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 4961a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 4971a989b97SToby Isaac 4981a989b97SToby Isaac if (j < Nminuskminus) { 499fad4db65SToby Isaac perm[l++] = i; 5001a989b97SToby Isaac Nk = Nminuskminus; 5011a989b97SToby Isaac } else { 5021a989b97SToby Isaac subcomp[m++] = i; 5031a989b97SToby Isaac j -= Nminuskminus; 5041a989b97SToby Isaac odd ^= ((k - l) & 1); 5051a989b97SToby Isaac Nk = Nminusk; 5061a989b97SToby Isaac } 5071a989b97SToby Isaac } 5085f80ce2aSJacob Faibussowitsch for (; i < n; i++) subcomp[m++] = i; 5091a989b97SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 5101a989b97SToby Isaac PetscFunctionReturn(0); 5111a989b97SToby Isaac } 5121a989b97SToby Isaac 513ef0bb6c7SMatthew G. Knepley struct _p_PetscTabulation { 514a5b23f4aSJose E. Roman PetscInt K; /* Indicates a k-jet, namely tabulated derivatives up to order k */ 51519815104SMartin Diehl PetscInt Nr; /* The number of tabulation replicas (often 1) */ 516ef0bb6c7SMatthew G. Knepley PetscInt Np; /* The number of tabulation points in a replica */ 517ef0bb6c7SMatthew G. Knepley PetscInt Nb; /* The number of functions tabulated */ 518ef0bb6c7SMatthew G. Knepley PetscInt Nc; /* The number of function components */ 519ef0bb6c7SMatthew G. Knepley PetscInt cdim; /* The coordinate dimension */ 520ef0bb6c7SMatthew G. Knepley PetscReal **T; /* The tabulation T[K] of functions and their derivatives 521ef0bb6c7SMatthew G. Knepley T[0] = B[Nr*Np][Nb][Nc]: The basis function values at quadrature points 522ef0bb6c7SMatthew G. Knepley T[1] = D[Nr*Np][Nb][Nc][cdim]: The basis function derivatives at quadrature points 523ef0bb6c7SMatthew G. Knepley T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */ 524ef0bb6c7SMatthew G. Knepley }; 525ef0bb6c7SMatthew G. Knepley typedef struct _p_PetscTabulation *PetscTabulation; 526ef0bb6c7SMatthew G. Knepley 527d6685f55SMatthew G. Knepley typedef PetscErrorCode (*PetscProbFunc)(const PetscReal[], const PetscReal[], PetscReal[]); 528d6685f55SMatthew G. Knepley 5299371c9d4SSatish Balay typedef enum { 5309371c9d4SSatish Balay DTPROB_DENSITY_CONSTANT, 5319371c9d4SSatish Balay DTPROB_DENSITY_GAUSSIAN, 5329371c9d4SSatish Balay DTPROB_DENSITY_MAXWELL_BOLTZMANN, 5339371c9d4SSatish Balay DTPROB_NUM_DENSITY 5349371c9d4SSatish Balay } DTProbDensityType; 535d6685f55SMatthew G. Knepley PETSC_EXTERN const char *const DTProbDensityTypes[]; 536d6685f55SMatthew G. Knepley 537d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]); 538d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]); 539d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]); 540d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]); 541d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]); 542d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]); 543d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 544d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 545d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 546d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]); 547d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]); 548ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]); 549ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]); 550d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 551d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 552d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 553ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 554ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 555ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 556ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 557ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 558ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 559d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFunc *, PetscProbFunc *, PetscProbFunc *); 560d6685f55SMatthew G. Knepley 561d6685f55SMatthew G. Knepley #include <petscvec.h> 562d6685f55SMatthew G. Knepley 563d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFunc, PetscReal *); 564d6685f55SMatthew G. Knepley 56537045ce4SJed Brown #endif 566