137045ce4SJed Brown /* 237045ce4SJed Brown Common tools for constructing discretizations 337045ce4SJed Brown */ 4*a4963045SJacob Faibussowitsch #pragma once 537045ce4SJed Brown 637045ce4SJed Brown #include <petscsys.h> 74366bac7SMatthew G. Knepley #include <petscdmtypes.h> 84366bac7SMatthew G. Knepley #include <petscistypes.h> 937045ce4SJed Brown 10ac09b921SBarry Smith /* SUBMANSEC = DT */ 11ac09b921SBarry Smith 122cd22861SMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCQUADRATURE_CLASSID; 132cd22861SMatthew G. Knepley 1421454ff5SMatthew G. Knepley /*S 1516a05f60SBarry Smith PetscQuadrature - Quadrature rule for numerical integration. 1621454ff5SMatthew G. Knepley 17329bbf4eSMatthew G. Knepley Level: beginner 1821454ff5SMatthew G. Knepley 19db781477SPatrick Sanan .seealso: `PetscQuadratureCreate()`, `PetscQuadratureDestroy()` 2021454ff5SMatthew G. Knepley S*/ 2121454ff5SMatthew G. Knepley typedef struct _p_PetscQuadrature *PetscQuadrature; 2221454ff5SMatthew G. Knepley 238272889dSSatish Balay /*E 24916e780bShannah_mairs PetscGaussLobattoLegendreCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights 258272889dSSatish Balay 2616a05f60SBarry Smith Values: 2716a05f60SBarry Smith + `PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA` - compute the nodes via linear algebra 2816a05f60SBarry Smith - `PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON` - compute the nodes by solving a nonlinear equation with Newton's method 2916a05f60SBarry Smith 308272889dSSatish Balay Level: intermediate 318272889dSSatish Balay 3216a05f60SBarry Smith .seealso: `PetscQuadrature` 338272889dSSatish Balay E*/ 349371c9d4SSatish Balay typedef enum { 359371c9d4SSatish Balay PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA, 369371c9d4SSatish Balay PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON 379371c9d4SSatish Balay } PetscGaussLobattoLegendreCreateType; 388272889dSSatish Balay 39d4afb720SToby Isaac /*E 40d4afb720SToby Isaac PetscDTNodeType - A description of strategies for generating nodes (both 41d4afb720SToby Isaac quadrature nodes and nodes for Lagrange polynomials) 42d4afb720SToby Isaac 4316a05f60SBarry Smith Values: 4416a05f60SBarry Smith + `PETSCDTNODES_DEFAULT` - Nodes chosen by PETSc 4516a05f60SBarry Smith . `PETSCDTNODES_GAUSSJACOBI` - Nodes at either Gauss-Jacobi or Gauss-Lobatto-Jacobi quadrature points 4616a05f60SBarry Smith . `PETSCDTNODES_EQUISPACED` - Nodes equispaced either including the endpoints or excluding them 4716a05f60SBarry Smith - `PETSCDTNODES_TANHSINH` - Nodes at Tanh-Sinh quadrature points 48d4afb720SToby Isaac 4916a05f60SBarry Smith Level: intermediate 50d4afb720SToby Isaac 5187497f52SBarry Smith Note: 5287497f52SBarry Smith A `PetscDTNodeType` can be paired with a `PetscBool` to indicate whether 5387497f52SBarry Smith the nodes include endpoints or not, and in the case of `PETSCDT_GAUSSJACOBI` 54d4afb720SToby Isaac with exponents for the weight function. 55d4afb720SToby Isaac 5616a05f60SBarry Smith .seealso: `PetscQuadrature` 57d4afb720SToby Isaac E*/ 589371c9d4SSatish Balay typedef enum { 599371c9d4SSatish Balay PETSCDTNODES_DEFAULT = -1, 609371c9d4SSatish Balay PETSCDTNODES_GAUSSJACOBI, 619371c9d4SSatish Balay PETSCDTNODES_EQUISPACED, 629371c9d4SSatish Balay PETSCDTNODES_TANHSINH 639371c9d4SSatish Balay } PetscDTNodeType; 64d4afb720SToby Isaac 65d3c69ad0SToby Isaac PETSC_EXTERN const char *const *const PetscDTNodeTypes; 66d3c69ad0SToby Isaac 67d3c69ad0SToby Isaac /*E 68d3c69ad0SToby Isaac PetscDTSimplexQuadratureType - A description of classes of quadrature rules for simplices 69d3c69ad0SToby Isaac 7016a05f60SBarry Smith Values: 7116a05f60SBarry Smith + `PETSCDTSIMPLEXQUAD_DEFAULT` - Quadrature rule chosen by PETSc 7216a05f60SBarry Smith . `PETSCDTSIMPLEXQUAD_CONIC` - Quadrature rules constructed as 73d3c69ad0SToby Isaac conically-warped tensor products of 1D 74d3c69ad0SToby Isaac Gauss-Jacobi quadrature rules. These are 75d3c69ad0SToby Isaac explicitly computable in any dimension for any 76d3c69ad0SToby Isaac degree, and the tensor-product structure can be 77d3c69ad0SToby Isaac exploited by sum-factorization methods, but 78d3c69ad0SToby Isaac they are not efficient in terms of nodes per 79d3c69ad0SToby Isaac polynomial degree. 8016a05f60SBarry Smith - `PETSCDTSIMPLEXQUAD_MINSYM` - Quadrature rules that are fully symmetric 81d3c69ad0SToby Isaac (symmetries of the simplex preserve the nodes 82d3c69ad0SToby Isaac and weights) with minimal (or near minimal) 83d3c69ad0SToby Isaac number of nodes. In dimensions higher than 1 84d3c69ad0SToby Isaac these are not simple to compute, so lookup 85d3c69ad0SToby Isaac tables are used. 86d3c69ad0SToby Isaac 8716a05f60SBarry Smith Level: intermediate 8816a05f60SBarry Smith 8916a05f60SBarry Smith .seealso: `PetscQuadrature`, `PetscDTSimplexQuadrature()` 90d3c69ad0SToby Isaac E*/ 919371c9d4SSatish Balay typedef enum { 929371c9d4SSatish Balay PETSCDTSIMPLEXQUAD_DEFAULT = -1, 939371c9d4SSatish Balay PETSCDTSIMPLEXQUAD_CONIC = 0, 949371c9d4SSatish Balay PETSCDTSIMPLEXQUAD_MINSYM 959371c9d4SSatish Balay } PetscDTSimplexQuadratureType; 96d3c69ad0SToby Isaac 97d3c69ad0SToby Isaac PETSC_EXTERN const char *const *const PetscDTSimplexQuadratureTypes; 98d4afb720SToby Isaac 9921454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *); 100c9638911SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *); 1014366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetCellType(PetscQuadrature, DMPolytopeType *); 1024366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetCellType(PetscQuadrature, DMPolytopeType); 103bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt *); 104bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt); 105a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt *); 106a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt); 1074f9ab2b4SJed Brown PETSC_EXTERN PetscErrorCode PetscQuadratureEqual(PetscQuadrature, PetscQuadrature, PetscBool *); 108a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt *, PetscInt *, PetscInt *, const PetscReal *[], const PetscReal *[]); 109a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal[], const PetscReal[]); 11021454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer); 11121454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *); 112a0845e3aSMatthew G. Knepley 1132df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTensorQuadratureCreate(PetscQuadrature, PetscQuadrature, PetscQuadrature *); 11489710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *); 1154366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureComputePermutations(PetscQuadrature, PetscInt *, IS *[]); 11689710940SMatthew G. Knepley 117907761f8SToby Isaac PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *); 118907761f8SToby Isaac 11937045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *); 120fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiNorm(PetscReal, PetscReal, PetscInt, PetscReal *); 12194e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt, PetscReal, PetscReal, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *); 122fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEvalJet(PetscReal, PetscReal, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]); 123fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPKDEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]); 124d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedSize(PetscInt, PetscInt, PetscInt, PetscInt *); 125d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscInt, PetscReal[]); 12637045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt, PetscReal, PetscReal, PetscReal *, PetscReal *); 12794e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *); 12894e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *); 129916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt, PetscGaussLobattoLegendreCreateType, PetscReal *, PetscReal *); 130194825f6SJed Brown PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *); 131a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *); 132e6a796c3SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *); 133d3c69ad0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTSimplexQuadrature(PetscInt, PetscInt, PetscDTSimplexQuadratureType, PetscQuadrature *); 1344366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTCreateDefaultQuadrature(DMPolytopeType, PetscInt, PetscQuadrature *, PetscQuadrature *); 13537045ce4SJed Brown 136b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *); 137d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *); 138d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *); 139b3c0f97bSTom Klotz 140916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *); 141916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 142916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 143916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 144916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 145916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 146916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 147916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 148916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 149916e780bShannah_mairs 1501a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 1511a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 1521a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 1531a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *); 1541a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *); 1551a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 1561a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *); 157dda711d0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]); 1581a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 1591a989b97SToby Isaac 160d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt, PetscInt, const PetscInt[], PetscInt *); 161d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt, PetscInt, PetscInt, PetscInt[]); 162fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGradedOrderToIndex(PetscInt, const PetscInt[], PetscInt *); 163fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToGradedOrder(PetscInt, PetscInt, PetscInt[]); 164d4afb720SToby Isaac 165fad4db65SToby Isaac #if defined(PETSC_USE_64BIT_INDICES) 166fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 20 167fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX 61 168fad4db65SToby Isaac #else 169fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 12 170fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX 29 171fad4db65SToby Isaac #endif 172fad4db65SToby Isaac 173fad4db65SToby Isaac /*MC 174fad4db65SToby Isaac PetscDTFactorial - Approximate n! as a real number 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 18316a05f60SBarry Smith 18416a05f60SBarry Smith .seealso: `PetscDTFactorialInt()`, `PetscDTBinomialInt()`, `PetscDTBinomial()` 185fad4db65SToby Isaac M*/ 186d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial) 187d71ae5a4SJacob Faibussowitsch { 188fad4db65SToby Isaac PetscReal f = 1.0; 189fad4db65SToby Isaac 190fad4db65SToby Isaac PetscFunctionBegin; 191e2ab39ccSLisandro Dalcin *factorial = -1.0; 19263a3b9bcSJacob Faibussowitsch PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %" PetscInt_FMT, n); 1935f80ce2aSJacob Faibussowitsch for (PetscInt i = 1; i < n + 1; ++i) f *= (PetscReal)i; 194fad4db65SToby Isaac *factorial = f; 1953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 196fad4db65SToby Isaac } 197fad4db65SToby Isaac 198fad4db65SToby Isaac /*MC 199fad4db65SToby Isaac PetscDTFactorialInt - Compute n! as an integer 200fad4db65SToby Isaac 2014165533cSJose E. Roman Input Parameter: 202fad4db65SToby Isaac . n - a non-negative integer 203fad4db65SToby Isaac 2044165533cSJose E. Roman Output Parameter: 205fad4db65SToby Isaac . factorial - n! 206fad4db65SToby Isaac 207fad4db65SToby Isaac Level: beginner 208fad4db65SToby Isaac 20916a05f60SBarry Smith Note: 21016a05f60SBarry Smith This is limited to n such that n! can be represented by `PetscInt`, which is 12 if `PetscInt` is a signed 32-bit integer and 20 if `PetscInt` is a signed 64-bit integer. 21116a05f60SBarry Smith 21216a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTBinomialInt()`, `PetscDTBinomial()` 213fad4db65SToby Isaac M*/ 214d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial) 215d71ae5a4SJacob Faibussowitsch { 216fad4db65SToby Isaac PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600}; 217fad4db65SToby Isaac 21828222859SToby Isaac PetscFunctionBegin; 21928222859SToby Isaac *factorial = -1; 22063a3b9bcSJacob 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); 221fad4db65SToby Isaac if (n <= 12) { 222fad4db65SToby Isaac *factorial = facLookup[n]; 223fad4db65SToby Isaac } else { 224fad4db65SToby Isaac PetscInt f = facLookup[12]; 225fad4db65SToby Isaac PetscInt i; 226fad4db65SToby Isaac 227fad4db65SToby Isaac for (i = 13; i < n + 1; ++i) f *= i; 228fad4db65SToby Isaac *factorial = f; 229fad4db65SToby Isaac } 2303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 231fad4db65SToby Isaac } 232fad4db65SToby Isaac 233fad4db65SToby Isaac /*MC 234fad4db65SToby Isaac PetscDTBinomial - Approximate the binomial coefficient "n choose k" 235fad4db65SToby Isaac 2364165533cSJose E. Roman Input Parameters: 237fad4db65SToby Isaac + n - a non-negative integer 238fad4db65SToby Isaac - k - an integer between 0 and n, inclusive 239fad4db65SToby Isaac 2404165533cSJose E. Roman Output Parameter: 241fad4db65SToby Isaac . binomial - approximation of the binomial coefficient n choose k 242fad4db65SToby Isaac 243fad4db65SToby Isaac Level: beginner 24416a05f60SBarry Smith 24516a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()` 246fad4db65SToby Isaac M*/ 247d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial) 248d71ae5a4SJacob Faibussowitsch { 2491a989b97SToby Isaac PetscFunctionBeginHot; 250e2ab39ccSLisandro Dalcin *binomial = -1.0; 25163a3b9bcSJacob 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); 2521a989b97SToby Isaac if (n <= 3) { 2539371c9d4SSatish Balay PetscInt binomLookup[4][4] = { 2549371c9d4SSatish Balay {1, 0, 0, 0}, 2559371c9d4SSatish Balay {1, 1, 0, 0}, 2569371c9d4SSatish Balay {1, 2, 1, 0}, 2579371c9d4SSatish Balay {1, 3, 3, 1} 2589371c9d4SSatish Balay }; 2591a989b97SToby Isaac 260e2ab39ccSLisandro Dalcin *binomial = (PetscReal)binomLookup[n][k]; 2611a989b97SToby Isaac } else { 262e2ab39ccSLisandro Dalcin PetscReal binom = 1.0; 2631a989b97SToby Isaac 2641a989b97SToby Isaac k = PetscMin(k, n - k); 2655f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1); 2661a989b97SToby Isaac *binomial = binom; 2671a989b97SToby Isaac } 2683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2691a989b97SToby Isaac } 2701a989b97SToby Isaac 271fad4db65SToby Isaac /*MC 272fad4db65SToby Isaac PetscDTBinomialInt - Compute the binomial coefficient "n choose k" 273fad4db65SToby Isaac 27497bb3fdcSJose E. Roman Input Parameters: 275fad4db65SToby Isaac + n - a non-negative integer 276fad4db65SToby Isaac - k - an integer between 0 and n, inclusive 277fad4db65SToby Isaac 27897bb3fdcSJose E. Roman Output Parameter: 279fad4db65SToby Isaac . binomial - the binomial coefficient n choose k 280fad4db65SToby Isaac 281fad4db65SToby Isaac Level: beginner 28216a05f60SBarry Smith 28316a05f60SBarry Smith Note: 28416a05f60SBarry Smith This is limited by integers that can be represented by `PetscInt`. 28516a05f60SBarry Smith 28616a05f60SBarry Smith Use `PetscDTBinomial()` for real number approximations of larger values 28716a05f60SBarry Smith 28816a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTEnumPerm()` 289fad4db65SToby Isaac M*/ 290d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial) 291d71ae5a4SJacob Faibussowitsch { 29228222859SToby Isaac PetscInt bin; 29328222859SToby Isaac 29428222859SToby Isaac PetscFunctionBegin; 29528222859SToby Isaac *binomial = -1; 29663a3b9bcSJacob 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); 29763a3b9bcSJacob 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); 298fad4db65SToby Isaac if (n <= 3) { 2999371c9d4SSatish Balay PetscInt binomLookup[4][4] = { 3009371c9d4SSatish Balay {1, 0, 0, 0}, 3019371c9d4SSatish Balay {1, 1, 0, 0}, 3029371c9d4SSatish Balay {1, 2, 1, 0}, 3039371c9d4SSatish Balay {1, 3, 3, 1} 3049371c9d4SSatish Balay }; 305fad4db65SToby Isaac 30628222859SToby Isaac bin = binomLookup[n][k]; 307fad4db65SToby Isaac } else { 308fad4db65SToby Isaac PetscInt binom = 1; 309fad4db65SToby Isaac 310fad4db65SToby Isaac k = PetscMin(k, n - k); 3115f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1); 31228222859SToby Isaac bin = binom; 313fad4db65SToby Isaac } 31428222859SToby Isaac *binomial = bin; 3153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 316fad4db65SToby Isaac } 317fad4db65SToby Isaac 318fad4db65SToby Isaac /*MC 31916a05f60SBarry Smith PetscDTEnumPerm - Get a permutation of `n` integers from its encoding into the integers [0, n!) as a sequence of swaps. 320fad4db65SToby Isaac 3214165533cSJose E. Roman Input Parameters: 322fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 3238cd1e013SToby Isaac - k - an integer in [0, n!) 324fad4db65SToby Isaac 3254165533cSJose E. Roman Output Parameters: 326fad4db65SToby Isaac + perm - the permuted list of the integers [0, ..., n-1] 3272fe279fdSBarry Smith - isOdd - if not `NULL`, returns whether the permutation used an even or odd number of swaps. 328fad4db65SToby Isaac 32916a05f60SBarry Smith Level: intermediate 330fad4db65SToby Isaac 33116a05f60SBarry Smith Notes: 33216a05f60SBarry Smith A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation, 33316a05f60SBarry Smith by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in 33416a05f60SBarry Smith some position j >= i. This swap is encoded as the difference (j - i). The difference d_i at step i is less than 33516a05f60SBarry Smith (n - i). This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number 33616a05f60SBarry Smith (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}. 33716a05f60SBarry Smith 33816a05f60SBarry 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. 33916a05f60SBarry Smith 34016a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTPermIndex()` 341fad4db65SToby Isaac M*/ 342d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd) 343d71ae5a4SJacob Faibussowitsch { 3441a989b97SToby Isaac PetscInt odd = 0; 3451a989b97SToby Isaac PetscInt i; 346fad4db65SToby Isaac PetscInt work[PETSC_FACTORIAL_MAX]; 347fad4db65SToby Isaac PetscInt *w; 3481a989b97SToby Isaac 34928222859SToby Isaac PetscFunctionBegin; 35028222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 35163a3b9bcSJacob 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); 352fad4db65SToby Isaac w = &work[n - 2]; 3531a989b97SToby Isaac for (i = 2; i <= n; i++) { 3541a989b97SToby Isaac *(w--) = k % i; 3551a989b97SToby Isaac k /= i; 3561a989b97SToby Isaac } 3571a989b97SToby Isaac for (i = 0; i < n; i++) perm[i] = i; 3581a989b97SToby Isaac for (i = 0; i < n - 1; i++) { 3591a989b97SToby Isaac PetscInt s = work[i]; 3601a989b97SToby Isaac PetscInt swap = perm[i]; 3611a989b97SToby Isaac 3621a989b97SToby Isaac perm[i] = perm[i + s]; 3631a989b97SToby Isaac perm[i + s] = swap; 3641a989b97SToby Isaac odd ^= (!!s); 3651a989b97SToby Isaac } 3661a989b97SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 3673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3681a989b97SToby Isaac } 3691a989b97SToby Isaac 370fad4db65SToby Isaac /*MC 37116a05f60SBarry Smith PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!). This inverts `PetscDTEnumPerm()`. 3728cd1e013SToby Isaac 3734165533cSJose E. Roman Input Parameters: 3748cd1e013SToby Isaac + n - a non-negative integer (see note about limits below) 3758cd1e013SToby Isaac - perm - the permuted list of the integers [0, ..., n-1] 3768cd1e013SToby Isaac 3774165533cSJose E. Roman Output Parameters: 3788cd1e013SToby Isaac + k - an integer in [0, n!) 3792fe279fdSBarry Smith - isOdd - if not `NULL`, returns whether the permutation used an even or odd number of swaps. 3808cd1e013SToby Isaac 3818cd1e013SToby Isaac Level: beginner 38216a05f60SBarry Smith 38316a05f60SBarry Smith Note: 38416a05f60SBarry 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. 38516a05f60SBarry Smith 38616a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()` 3878cd1e013SToby Isaac M*/ 388d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd) 389d71ae5a4SJacob Faibussowitsch { 3908cd1e013SToby Isaac PetscInt odd = 0; 3918cd1e013SToby Isaac PetscInt i, idx; 3928cd1e013SToby Isaac PetscInt work[PETSC_FACTORIAL_MAX]; 3938cd1e013SToby Isaac PetscInt iwork[PETSC_FACTORIAL_MAX]; 3948cd1e013SToby Isaac 3958cd1e013SToby Isaac PetscFunctionBeginHot; 39628222859SToby Isaac *k = -1; 39728222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 39863a3b9bcSJacob 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); 3998cd1e013SToby Isaac for (i = 0; i < n; i++) work[i] = i; /* partial permutation */ 4008cd1e013SToby Isaac for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */ 4018cd1e013SToby Isaac for (idx = 0, i = 0; i < n - 1; i++) { 4028cd1e013SToby Isaac PetscInt j = perm[i]; 4038cd1e013SToby Isaac PetscInt icur = work[i]; 4048cd1e013SToby Isaac PetscInt jloc = iwork[j]; 4058cd1e013SToby Isaac PetscInt diff = jloc - i; 4068cd1e013SToby Isaac 4078cd1e013SToby Isaac idx = idx * (n - i) + diff; 4088cd1e013SToby Isaac /* swap (i, jloc) */ 4098cd1e013SToby Isaac work[i] = j; 4108cd1e013SToby Isaac work[jloc] = icur; 4118cd1e013SToby Isaac iwork[j] = i; 4128cd1e013SToby Isaac iwork[icur] = jloc; 4138cd1e013SToby Isaac odd ^= (!!diff); 4148cd1e013SToby Isaac } 4158cd1e013SToby Isaac *k = idx; 4168cd1e013SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 4173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4188cd1e013SToby Isaac } 4198cd1e013SToby Isaac 4208cd1e013SToby Isaac /*MC 421fad4db65SToby Isaac PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k). 422fad4db65SToby Isaac The encoding is in lexicographic order. 423fad4db65SToby Isaac 4244165533cSJose E. Roman Input Parameters: 425fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 426fad4db65SToby Isaac . k - an integer in [0, n] 427fad4db65SToby Isaac - j - an index in [0, n choose k) 428fad4db65SToby Isaac 4294165533cSJose E. Roman Output Parameter: 430fad4db65SToby Isaac . subset - the jth subset of size k of the integers [0, ..., n - 1] 431fad4db65SToby Isaac 432fad4db65SToby Isaac Level: beginner 433fad4db65SToby Isaac 43416a05f60SBarry Smith Note: 43516a05f60SBarry Smith Limited by arguments such that `n` choose `k` can be represented by `PetscInt` 43616a05f60SBarry Smith 43716a05f60SBarry Smith .seealso: `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()` 438fad4db65SToby Isaac M*/ 439d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset) 440d71ae5a4SJacob Faibussowitsch { 4415f80ce2aSJacob Faibussowitsch PetscInt Nk; 4421a989b97SToby Isaac 4431a989b97SToby Isaac PetscFunctionBeginHot; 4449566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 4455f80ce2aSJacob Faibussowitsch for (PetscInt i = 0, l = 0; i < n && l < k; i++) { 4461a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 4471a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 4481a989b97SToby Isaac 4491a989b97SToby Isaac if (j < Nminuskminus) { 4501a989b97SToby Isaac subset[l++] = i; 4511a989b97SToby Isaac Nk = Nminuskminus; 4521a989b97SToby Isaac } else { 4531a989b97SToby Isaac j -= Nminuskminus; 4541a989b97SToby Isaac Nk = Nminusk; 4551a989b97SToby Isaac } 4561a989b97SToby Isaac } 4573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4581a989b97SToby Isaac } 4591a989b97SToby Isaac 460fad4db65SToby Isaac /*MC 46187497f52SBarry 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. 46287497f52SBarry Smith This is the inverse of `PetscDTEnumSubset`. 463fad4db65SToby Isaac 4644165533cSJose E. Roman Input Parameters: 465fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 466fad4db65SToby Isaac . k - an integer in [0, n] 467fad4db65SToby Isaac - subset - an ordered subset of the integers [0, ..., n - 1] 468fad4db65SToby Isaac 4694165533cSJose E. Roman Output Parameter: 470fad4db65SToby Isaac . index - the rank of the subset in lexicographic order 471fad4db65SToby Isaac 472fad4db65SToby Isaac Level: beginner 473fad4db65SToby Isaac 47416a05f60SBarry Smith Note: 47516a05f60SBarry Smith Limited by arguments such that `n` choose `k` can be represented by `PetscInt` 47616a05f60SBarry Smith 47716a05f60SBarry Smith .seealso: `PetscDTEnumSubset()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()` 478fad4db65SToby Isaac M*/ 479d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index) 480d71ae5a4SJacob Faibussowitsch { 4815f80ce2aSJacob Faibussowitsch PetscInt j = 0, Nk; 4821a989b97SToby Isaac 48328222859SToby Isaac PetscFunctionBegin; 48428222859SToby Isaac *index = -1; 4859566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 4865f80ce2aSJacob Faibussowitsch for (PetscInt i = 0, l = 0; i < n && l < k; i++) { 4871a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 4881a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 4891a989b97SToby Isaac 4901a989b97SToby Isaac if (subset[l] == i) { 4911a989b97SToby Isaac l++; 4921a989b97SToby Isaac Nk = Nminuskminus; 4931a989b97SToby Isaac } else { 4941a989b97SToby Isaac j += Nminuskminus; 4951a989b97SToby Isaac Nk = Nminusk; 4961a989b97SToby Isaac } 4971a989b97SToby Isaac } 4981a989b97SToby Isaac *index = j; 4993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5001a989b97SToby Isaac } 5011a989b97SToby Isaac 502fad4db65SToby Isaac /*MC 50328222859SToby 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. 504fad4db65SToby Isaac 5054165533cSJose E. Roman Input Parameters: 506fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 507fad4db65SToby Isaac . k - an integer in [0, n] 508fad4db65SToby Isaac - j - an index in [0, n choose k) 509fad4db65SToby Isaac 5104165533cSJose E. Roman Output Parameters: 511fad4db65SToby Isaac + perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set. 5122fe279fdSBarry Smith - isOdd - if not `NULL`, return whether perm is an even or odd permutation. 513fad4db65SToby Isaac 514fad4db65SToby Isaac Level: beginner 515fad4db65SToby Isaac 51616a05f60SBarry Smith Note: 51716a05f60SBarry Smith Limited by arguments such that `n` choose `k` can be represented by `PetscInt` 51816a05f60SBarry Smith 51916a05f60SBarry Smith .seealso: `PetscDTEnumSubset()`, `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, 52016a05f60SBarry Smith `PetscDTPermIndex()` 521fad4db65SToby Isaac M*/ 522d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd) 523d71ae5a4SJacob Faibussowitsch { 5245f80ce2aSJacob Faibussowitsch PetscInt i, l, m, Nk, odd = 0; 5255f80ce2aSJacob Faibussowitsch PetscInt *subcomp = perm + k; 5261a989b97SToby Isaac 52728222859SToby Isaac PetscFunctionBegin; 52828222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 5299566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 5301a989b97SToby Isaac for (i = 0, l = 0, m = 0; i < n && l < k; i++) { 5311a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 5321a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 5331a989b97SToby Isaac 5341a989b97SToby Isaac if (j < Nminuskminus) { 535fad4db65SToby Isaac perm[l++] = i; 5361a989b97SToby Isaac Nk = Nminuskminus; 5371a989b97SToby Isaac } else { 5381a989b97SToby Isaac subcomp[m++] = i; 5391a989b97SToby Isaac j -= Nminuskminus; 5401a989b97SToby Isaac odd ^= ((k - l) & 1); 5411a989b97SToby Isaac Nk = Nminusk; 5421a989b97SToby Isaac } 5431a989b97SToby Isaac } 5445f80ce2aSJacob Faibussowitsch for (; i < n; i++) subcomp[m++] = i; 5451a989b97SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 5463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5471a989b97SToby Isaac } 5481a989b97SToby Isaac 549ef0bb6c7SMatthew G. Knepley struct _p_PetscTabulation { 550a5b23f4aSJose E. Roman PetscInt K; /* Indicates a k-jet, namely tabulated derivatives up to order k */ 55119815104SMartin Diehl PetscInt Nr; /* The number of tabulation replicas (often 1) */ 552ef0bb6c7SMatthew G. Knepley PetscInt Np; /* The number of tabulation points in a replica */ 553ef0bb6c7SMatthew G. Knepley PetscInt Nb; /* The number of functions tabulated */ 554ef0bb6c7SMatthew G. Knepley PetscInt Nc; /* The number of function components */ 555ef0bb6c7SMatthew G. Knepley PetscInt cdim; /* The coordinate dimension */ 556ef0bb6c7SMatthew G. Knepley PetscReal **T; /* The tabulation T[K] of functions and their derivatives 557ef0bb6c7SMatthew G. Knepley T[0] = B[Nr*Np][Nb][Nc]: The basis function values at quadrature points 558ef0bb6c7SMatthew G. Knepley T[1] = D[Nr*Np][Nb][Nc][cdim]: The basis function derivatives at quadrature points 559ef0bb6c7SMatthew G. Knepley T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */ 560ef0bb6c7SMatthew G. Knepley }; 561ef0bb6c7SMatthew G. Knepley typedef struct _p_PetscTabulation *PetscTabulation; 562ef0bb6c7SMatthew G. Knepley 563d6685f55SMatthew G. Knepley typedef PetscErrorCode (*PetscProbFunc)(const PetscReal[], const PetscReal[], PetscReal[]); 564d6685f55SMatthew G. Knepley 5659371c9d4SSatish Balay typedef enum { 5669371c9d4SSatish Balay DTPROB_DENSITY_CONSTANT, 5679371c9d4SSatish Balay DTPROB_DENSITY_GAUSSIAN, 5689371c9d4SSatish Balay DTPROB_DENSITY_MAXWELL_BOLTZMANN, 5699371c9d4SSatish Balay DTPROB_NUM_DENSITY 5709371c9d4SSatish Balay } DTProbDensityType; 571d6685f55SMatthew G. Knepley PETSC_EXTERN const char *const DTProbDensityTypes[]; 572d6685f55SMatthew G. Knepley 573d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]); 574d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]); 575d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]); 576d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]); 577d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]); 578d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]); 579d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 580d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 581d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 582d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]); 583d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]); 584ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]); 585ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]); 586d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 587d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 588d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 589ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 590ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 591ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 592ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 593ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 594ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 595d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFunc *, PetscProbFunc *, PetscProbFunc *); 596d6685f55SMatthew G. Knepley 597d6685f55SMatthew G. Knepley #include <petscvec.h> 598d6685f55SMatthew G. Knepley 599d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFunc, PetscReal *); 600