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> 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*/ 174d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial) 175d71ae5a4SJacob Faibussowitsch { 176fad4db65SToby Isaac PetscReal f = 1.0; 177fad4db65SToby Isaac 178fad4db65SToby Isaac PetscFunctionBegin; 179e2ab39ccSLisandro Dalcin *factorial = -1.0; 18063a3b9bcSJacob Faibussowitsch PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %" PetscInt_FMT, n); 1815f80ce2aSJacob Faibussowitsch for (PetscInt i = 1; i < n + 1; ++i) f *= (PetscReal)i; 182fad4db65SToby Isaac *factorial = f; 1833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 184fad4db65SToby Isaac } 185fad4db65SToby Isaac 186fad4db65SToby Isaac /*MC 187fad4db65SToby Isaac PetscDTFactorialInt - Compute n! as an integer 188fad4db65SToby Isaac 1894165533cSJose E. Roman Input Parameter: 190fad4db65SToby Isaac . n - a non-negative integer 191fad4db65SToby Isaac 1924165533cSJose E. Roman Output Parameter: 193fad4db65SToby Isaac . factorial - n! 194fad4db65SToby Isaac 195fad4db65SToby Isaac Level: beginner 196fad4db65SToby Isaac 197fad4db65SToby 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. 198fad4db65SToby Isaac M*/ 199d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial) 200d71ae5a4SJacob Faibussowitsch { 201fad4db65SToby Isaac PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600}; 202fad4db65SToby Isaac 20328222859SToby Isaac PetscFunctionBegin; 20428222859SToby Isaac *factorial = -1; 20563a3b9bcSJacob 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); 206fad4db65SToby Isaac if (n <= 12) { 207fad4db65SToby Isaac *factorial = facLookup[n]; 208fad4db65SToby Isaac } else { 209fad4db65SToby Isaac PetscInt f = facLookup[12]; 210fad4db65SToby Isaac PetscInt i; 211fad4db65SToby Isaac 212fad4db65SToby Isaac for (i = 13; i < n + 1; ++i) f *= i; 213fad4db65SToby Isaac *factorial = f; 214fad4db65SToby Isaac } 2153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 216fad4db65SToby Isaac } 217fad4db65SToby Isaac 218fad4db65SToby Isaac /*MC 219fad4db65SToby Isaac PetscDTBinomial - Approximate the binomial coefficient "n choose k" 220fad4db65SToby Isaac 2214165533cSJose E. Roman Input Parameters: 222fad4db65SToby Isaac + n - a non-negative integer 223fad4db65SToby Isaac - k - an integer between 0 and n, inclusive 224fad4db65SToby Isaac 2254165533cSJose E. Roman Output Parameter: 226fad4db65SToby Isaac . binomial - approximation of the binomial coefficient n choose k 227fad4db65SToby Isaac 228fad4db65SToby Isaac Level: beginner 229fad4db65SToby Isaac M*/ 230d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial) 231d71ae5a4SJacob Faibussowitsch { 2321a989b97SToby Isaac PetscFunctionBeginHot; 233e2ab39ccSLisandro Dalcin *binomial = -1.0; 23463a3b9bcSJacob 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); 2351a989b97SToby Isaac if (n <= 3) { 2369371c9d4SSatish Balay PetscInt binomLookup[4][4] = { 2379371c9d4SSatish Balay {1, 0, 0, 0}, 2389371c9d4SSatish Balay {1, 1, 0, 0}, 2399371c9d4SSatish Balay {1, 2, 1, 0}, 2409371c9d4SSatish Balay {1, 3, 3, 1} 2419371c9d4SSatish Balay }; 2421a989b97SToby Isaac 243e2ab39ccSLisandro Dalcin *binomial = (PetscReal)binomLookup[n][k]; 2441a989b97SToby Isaac } else { 245e2ab39ccSLisandro Dalcin PetscReal binom = 1.0; 2461a989b97SToby Isaac 2471a989b97SToby Isaac k = PetscMin(k, n - k); 2485f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1); 2491a989b97SToby Isaac *binomial = binom; 2501a989b97SToby Isaac } 2513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2521a989b97SToby Isaac } 2531a989b97SToby Isaac 254fad4db65SToby Isaac /*MC 255fad4db65SToby Isaac PetscDTBinomialInt - Compute the binomial coefficient "n choose k" 256fad4db65SToby Isaac 25797bb3fdcSJose E. Roman Input Parameters: 258fad4db65SToby Isaac + n - a non-negative integer 259fad4db65SToby Isaac - k - an integer between 0 and n, inclusive 260fad4db65SToby Isaac 26197bb3fdcSJose E. Roman Output Parameter: 262fad4db65SToby Isaac . binomial - the binomial coefficient n choose k 263fad4db65SToby Isaac 26487497f52SBarry Smith Note: this is limited by integers that can be represented by `PetscInt` 265fad4db65SToby Isaac 266fad4db65SToby Isaac Level: beginner 267fad4db65SToby Isaac M*/ 268d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial) 269d71ae5a4SJacob Faibussowitsch { 27028222859SToby Isaac PetscInt bin; 27128222859SToby Isaac 27228222859SToby Isaac PetscFunctionBegin; 27328222859SToby Isaac *binomial = -1; 27463a3b9bcSJacob 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); 27563a3b9bcSJacob 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); 276fad4db65SToby Isaac if (n <= 3) { 2779371c9d4SSatish Balay PetscInt binomLookup[4][4] = { 2789371c9d4SSatish Balay {1, 0, 0, 0}, 2799371c9d4SSatish Balay {1, 1, 0, 0}, 2809371c9d4SSatish Balay {1, 2, 1, 0}, 2819371c9d4SSatish Balay {1, 3, 3, 1} 2829371c9d4SSatish Balay }; 283fad4db65SToby Isaac 28428222859SToby Isaac bin = binomLookup[n][k]; 285fad4db65SToby Isaac } else { 286fad4db65SToby Isaac PetscInt binom = 1; 287fad4db65SToby Isaac 288fad4db65SToby Isaac k = PetscMin(k, n - k); 2895f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1); 29028222859SToby Isaac bin = binom; 291fad4db65SToby Isaac } 29228222859SToby Isaac *binomial = bin; 2933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 294fad4db65SToby Isaac } 295fad4db65SToby Isaac 296fad4db65SToby Isaac /*MC 297fad4db65SToby Isaac PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps. 298fad4db65SToby Isaac 299fad4db65SToby Isaac A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation, 300fad4db65SToby Isaac by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in 30128222859SToby Isaac some position j >= i. This swap is encoded as the difference (j - i). The difference d_i at step i is less than 30228222859SToby Isaac (n - i). This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number 303fad4db65SToby Isaac (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}. 304fad4db65SToby Isaac 3054165533cSJose E. Roman Input Parameters: 306fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 3078cd1e013SToby Isaac - k - an integer in [0, n!) 308fad4db65SToby Isaac 3094165533cSJose E. Roman Output Parameters: 310fad4db65SToby Isaac + perm - the permuted list of the integers [0, ..., n-1] 311*2fe279fdSBarry Smith - isOdd - if not `NULL`, returns whether the permutation used an even or odd number of swaps. 312fad4db65SToby Isaac 31387497f52SBarry Smith Note: 31487497f52SBarry 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. 315fad4db65SToby Isaac 316fad4db65SToby Isaac Level: beginner 317fad4db65SToby Isaac M*/ 318d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd) 319d71ae5a4SJacob Faibussowitsch { 3201a989b97SToby Isaac PetscInt odd = 0; 3211a989b97SToby Isaac PetscInt i; 322fad4db65SToby Isaac PetscInt work[PETSC_FACTORIAL_MAX]; 323fad4db65SToby Isaac PetscInt *w; 3241a989b97SToby Isaac 32528222859SToby Isaac PetscFunctionBegin; 32628222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 32763a3b9bcSJacob 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); 328fad4db65SToby Isaac w = &work[n - 2]; 3291a989b97SToby Isaac for (i = 2; i <= n; i++) { 3301a989b97SToby Isaac *(w--) = k % i; 3311a989b97SToby Isaac k /= i; 3321a989b97SToby Isaac } 3331a989b97SToby Isaac for (i = 0; i < n; i++) perm[i] = i; 3341a989b97SToby Isaac for (i = 0; i < n - 1; i++) { 3351a989b97SToby Isaac PetscInt s = work[i]; 3361a989b97SToby Isaac PetscInt swap = perm[i]; 3371a989b97SToby Isaac 3381a989b97SToby Isaac perm[i] = perm[i + s]; 3391a989b97SToby Isaac perm[i + s] = swap; 3401a989b97SToby Isaac odd ^= (!!s); 3411a989b97SToby Isaac } 3421a989b97SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 3433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3441a989b97SToby Isaac } 3451a989b97SToby Isaac 346fad4db65SToby Isaac /*MC 34787497f52SBarry Smith PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!). This inverts `PetscDTEnumPerm`. 3488cd1e013SToby Isaac 3494165533cSJose E. Roman Input Parameters: 3508cd1e013SToby Isaac + n - a non-negative integer (see note about limits below) 3518cd1e013SToby Isaac - perm - the permuted list of the integers [0, ..., n-1] 3528cd1e013SToby Isaac 3534165533cSJose E. Roman Output Parameters: 3548cd1e013SToby Isaac + k - an integer in [0, n!) 355*2fe279fdSBarry Smith - isOdd - if not `NULL`, returns whether the permutation used an even or odd number of swaps. 3568cd1e013SToby Isaac 35787497f52SBarry Smith Note: 35887497f52SBarry 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. 3598cd1e013SToby Isaac 3608cd1e013SToby Isaac Level: beginner 3618cd1e013SToby Isaac M*/ 362d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd) 363d71ae5a4SJacob Faibussowitsch { 3648cd1e013SToby Isaac PetscInt odd = 0; 3658cd1e013SToby Isaac PetscInt i, idx; 3668cd1e013SToby Isaac PetscInt work[PETSC_FACTORIAL_MAX]; 3678cd1e013SToby Isaac PetscInt iwork[PETSC_FACTORIAL_MAX]; 3688cd1e013SToby Isaac 3698cd1e013SToby Isaac PetscFunctionBeginHot; 37028222859SToby Isaac *k = -1; 37128222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 37263a3b9bcSJacob 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); 3738cd1e013SToby Isaac for (i = 0; i < n; i++) work[i] = i; /* partial permutation */ 3748cd1e013SToby Isaac for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */ 3758cd1e013SToby Isaac for (idx = 0, i = 0; i < n - 1; i++) { 3768cd1e013SToby Isaac PetscInt j = perm[i]; 3778cd1e013SToby Isaac PetscInt icur = work[i]; 3788cd1e013SToby Isaac PetscInt jloc = iwork[j]; 3798cd1e013SToby Isaac PetscInt diff = jloc - i; 3808cd1e013SToby Isaac 3818cd1e013SToby Isaac idx = idx * (n - i) + diff; 3828cd1e013SToby Isaac /* swap (i, jloc) */ 3838cd1e013SToby Isaac work[i] = j; 3848cd1e013SToby Isaac work[jloc] = icur; 3858cd1e013SToby Isaac iwork[j] = i; 3868cd1e013SToby Isaac iwork[icur] = jloc; 3878cd1e013SToby Isaac odd ^= (!!diff); 3888cd1e013SToby Isaac } 3898cd1e013SToby Isaac *k = idx; 3908cd1e013SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 3913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3928cd1e013SToby Isaac } 3938cd1e013SToby Isaac 3948cd1e013SToby Isaac /*MC 395fad4db65SToby Isaac PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k). 396fad4db65SToby Isaac The encoding is in lexicographic order. 397fad4db65SToby Isaac 3984165533cSJose E. Roman Input Parameters: 399fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 400fad4db65SToby Isaac . k - an integer in [0, n] 401fad4db65SToby Isaac - j - an index in [0, n choose k) 402fad4db65SToby Isaac 4034165533cSJose E. Roman Output Parameter: 404fad4db65SToby Isaac . subset - the jth subset of size k of the integers [0, ..., n - 1] 405fad4db65SToby Isaac 40687497f52SBarry Smith Note: 40787497f52SBarry Smith Limited by arguments such that n choose k can be represented by `PetscInt` 408fad4db65SToby Isaac 409fad4db65SToby Isaac Level: beginner 410fad4db65SToby Isaac 411db781477SPatrick Sanan .seealso: `PetscDTSubsetIndex()` 412fad4db65SToby Isaac M*/ 413d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset) 414d71ae5a4SJacob Faibussowitsch { 4155f80ce2aSJacob Faibussowitsch PetscInt Nk; 4161a989b97SToby Isaac 4171a989b97SToby Isaac PetscFunctionBeginHot; 4189566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 4195f80ce2aSJacob Faibussowitsch for (PetscInt i = 0, l = 0; i < n && l < k; i++) { 4201a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 4211a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 4221a989b97SToby Isaac 4231a989b97SToby Isaac if (j < Nminuskminus) { 4241a989b97SToby Isaac subset[l++] = i; 4251a989b97SToby Isaac Nk = Nminuskminus; 4261a989b97SToby Isaac } else { 4271a989b97SToby Isaac j -= Nminuskminus; 4281a989b97SToby Isaac Nk = Nminusk; 4291a989b97SToby Isaac } 4301a989b97SToby Isaac } 4313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4321a989b97SToby Isaac } 4331a989b97SToby Isaac 434fad4db65SToby Isaac /*MC 43587497f52SBarry 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. 43687497f52SBarry Smith This is the inverse of `PetscDTEnumSubset`. 437fad4db65SToby Isaac 4384165533cSJose E. Roman Input Parameters: 439fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 440fad4db65SToby Isaac . k - an integer in [0, n] 441fad4db65SToby Isaac - subset - an ordered subset of the integers [0, ..., n - 1] 442fad4db65SToby Isaac 4434165533cSJose E. Roman Output Parameter: 444fad4db65SToby Isaac . index - the rank of the subset in lexicographic order 445fad4db65SToby Isaac 44687497f52SBarry Smith Note: 44787497f52SBarry Smith Limited by arguments such that n choose k can be represented by `PetscInt` 448fad4db65SToby Isaac 449fad4db65SToby Isaac Level: beginner 450fad4db65SToby Isaac 451db781477SPatrick Sanan .seealso: `PetscDTEnumSubset()` 452fad4db65SToby Isaac M*/ 453d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index) 454d71ae5a4SJacob Faibussowitsch { 4555f80ce2aSJacob Faibussowitsch PetscInt j = 0, Nk; 4561a989b97SToby Isaac 45728222859SToby Isaac PetscFunctionBegin; 45828222859SToby Isaac *index = -1; 4599566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 4605f80ce2aSJacob Faibussowitsch for (PetscInt i = 0, l = 0; i < n && l < k; i++) { 4611a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 4621a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 4631a989b97SToby Isaac 4641a989b97SToby Isaac if (subset[l] == i) { 4651a989b97SToby Isaac l++; 4661a989b97SToby Isaac Nk = Nminuskminus; 4671a989b97SToby Isaac } else { 4681a989b97SToby Isaac j += Nminuskminus; 4691a989b97SToby Isaac Nk = Nminusk; 4701a989b97SToby Isaac } 4711a989b97SToby Isaac } 4721a989b97SToby Isaac *index = j; 4733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4741a989b97SToby Isaac } 4751a989b97SToby Isaac 476fad4db65SToby Isaac /*MC 47728222859SToby 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. 478fad4db65SToby Isaac 4794165533cSJose E. Roman Input Parameters: 480fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 481fad4db65SToby Isaac . k - an integer in [0, n] 482fad4db65SToby Isaac - j - an index in [0, n choose k) 483fad4db65SToby Isaac 4844165533cSJose E. Roman Output Parameters: 485fad4db65SToby Isaac + perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set. 486*2fe279fdSBarry Smith - isOdd - if not `NULL`, return whether perm is an even or odd permutation. 487fad4db65SToby Isaac 48887497f52SBarry Smith Note: 48987497f52SBarry Smith Limited by arguments such that n choose k can be represented by `PetscInt` 490fad4db65SToby Isaac 491fad4db65SToby Isaac Level: beginner 492fad4db65SToby Isaac 493db781477SPatrick Sanan .seealso: `PetscDTEnumSubset()`, `PetscDTSubsetIndex()` 494fad4db65SToby Isaac M*/ 495d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd) 496d71ae5a4SJacob Faibussowitsch { 4975f80ce2aSJacob Faibussowitsch PetscInt i, l, m, Nk, odd = 0; 4985f80ce2aSJacob Faibussowitsch PetscInt *subcomp = perm + k; 4991a989b97SToby Isaac 50028222859SToby Isaac PetscFunctionBegin; 50128222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 5029566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 5031a989b97SToby Isaac for (i = 0, l = 0, m = 0; i < n && l < k; i++) { 5041a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 5051a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 5061a989b97SToby Isaac 5071a989b97SToby Isaac if (j < Nminuskminus) { 508fad4db65SToby Isaac perm[l++] = i; 5091a989b97SToby Isaac Nk = Nminuskminus; 5101a989b97SToby Isaac } else { 5111a989b97SToby Isaac subcomp[m++] = i; 5121a989b97SToby Isaac j -= Nminuskminus; 5131a989b97SToby Isaac odd ^= ((k - l) & 1); 5141a989b97SToby Isaac Nk = Nminusk; 5151a989b97SToby Isaac } 5161a989b97SToby Isaac } 5175f80ce2aSJacob Faibussowitsch for (; i < n; i++) subcomp[m++] = i; 5181a989b97SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 5193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5201a989b97SToby Isaac } 5211a989b97SToby Isaac 522ef0bb6c7SMatthew G. Knepley struct _p_PetscTabulation { 523a5b23f4aSJose E. Roman PetscInt K; /* Indicates a k-jet, namely tabulated derivatives up to order k */ 52419815104SMartin Diehl PetscInt Nr; /* The number of tabulation replicas (often 1) */ 525ef0bb6c7SMatthew G. Knepley PetscInt Np; /* The number of tabulation points in a replica */ 526ef0bb6c7SMatthew G. Knepley PetscInt Nb; /* The number of functions tabulated */ 527ef0bb6c7SMatthew G. Knepley PetscInt Nc; /* The number of function components */ 528ef0bb6c7SMatthew G. Knepley PetscInt cdim; /* The coordinate dimension */ 529ef0bb6c7SMatthew G. Knepley PetscReal **T; /* The tabulation T[K] of functions and their derivatives 530ef0bb6c7SMatthew G. Knepley T[0] = B[Nr*Np][Nb][Nc]: The basis function values at quadrature points 531ef0bb6c7SMatthew G. Knepley T[1] = D[Nr*Np][Nb][Nc][cdim]: The basis function derivatives at quadrature points 532ef0bb6c7SMatthew G. Knepley T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */ 533ef0bb6c7SMatthew G. Knepley }; 534ef0bb6c7SMatthew G. Knepley typedef struct _p_PetscTabulation *PetscTabulation; 535ef0bb6c7SMatthew G. Knepley 536d6685f55SMatthew G. Knepley typedef PetscErrorCode (*PetscProbFunc)(const PetscReal[], const PetscReal[], PetscReal[]); 537d6685f55SMatthew G. Knepley 5389371c9d4SSatish Balay typedef enum { 5399371c9d4SSatish Balay DTPROB_DENSITY_CONSTANT, 5409371c9d4SSatish Balay DTPROB_DENSITY_GAUSSIAN, 5419371c9d4SSatish Balay DTPROB_DENSITY_MAXWELL_BOLTZMANN, 5429371c9d4SSatish Balay DTPROB_NUM_DENSITY 5439371c9d4SSatish Balay } DTProbDensityType; 544d6685f55SMatthew G. Knepley PETSC_EXTERN const char *const DTProbDensityTypes[]; 545d6685f55SMatthew G. Knepley 546d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]); 547d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]); 548d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]); 549d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]); 550d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]); 551d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]); 552d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 553d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 554d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 555d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]); 556d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]); 557ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]); 558ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]); 559d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 560d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 561d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 562ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 563ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 564ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 565ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 566ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 567ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 568d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFunc *, PetscProbFunc *, PetscProbFunc *); 569d6685f55SMatthew G. Knepley 570d6685f55SMatthew G. Knepley #include <petscvec.h> 571d6685f55SMatthew G. Knepley 572d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFunc, PetscReal *); 573d6685f55SMatthew G. Knepley 57437045ce4SJed Brown #endif 575