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> 84366bac7SMatthew G. Knepley #include <petscdmtypes.h> 94366bac7SMatthew G. Knepley #include <petscistypes.h> 1037045ce4SJed Brown 11ac09b921SBarry Smith /* SUBMANSEC = DT */ 12ac09b921SBarry Smith 132cd22861SMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCQUADRATURE_CLASSID; 142cd22861SMatthew G. Knepley 1521454ff5SMatthew G. Knepley /*S 16*16a05f60SBarry Smith PetscQuadrature - Quadrature rule for numerical integration. 1721454ff5SMatthew G. Knepley 18329bbf4eSMatthew G. Knepley Level: beginner 1921454ff5SMatthew G. Knepley 20db781477SPatrick Sanan .seealso: `PetscQuadratureCreate()`, `PetscQuadratureDestroy()` 2121454ff5SMatthew G. Knepley S*/ 2221454ff5SMatthew G. Knepley typedef struct _p_PetscQuadrature *PetscQuadrature; 2321454ff5SMatthew G. Knepley 248272889dSSatish Balay /*E 25916e780bShannah_mairs PetscGaussLobattoLegendreCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights 268272889dSSatish Balay 27*16a05f60SBarry Smith Values: 28*16a05f60SBarry Smith + `PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA` - compute the nodes via linear algebra 29*16a05f60SBarry Smith - `PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON` - compute the nodes by solving a nonlinear equation with Newton's method 30*16a05f60SBarry Smith 318272889dSSatish Balay Level: intermediate 328272889dSSatish Balay 33*16a05f60SBarry Smith .seealso: `PetscQuadrature` 348272889dSSatish Balay E*/ 359371c9d4SSatish Balay typedef enum { 369371c9d4SSatish Balay PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA, 379371c9d4SSatish Balay PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON 389371c9d4SSatish Balay } PetscGaussLobattoLegendreCreateType; 398272889dSSatish Balay 40d4afb720SToby Isaac /*E 41d4afb720SToby Isaac PetscDTNodeType - A description of strategies for generating nodes (both 42d4afb720SToby Isaac quadrature nodes and nodes for Lagrange polynomials) 43d4afb720SToby Isaac 44*16a05f60SBarry Smith Values: 45*16a05f60SBarry Smith + `PETSCDTNODES_DEFAULT` - Nodes chosen by PETSc 46*16a05f60SBarry Smith . `PETSCDTNODES_GAUSSJACOBI` - Nodes at either Gauss-Jacobi or Gauss-Lobatto-Jacobi quadrature points 47*16a05f60SBarry Smith . `PETSCDTNODES_EQUISPACED` - Nodes equispaced either including the endpoints or excluding them 48*16a05f60SBarry Smith - `PETSCDTNODES_TANHSINH` - Nodes at Tanh-Sinh quadrature points 49d4afb720SToby Isaac 50*16a05f60SBarry Smith Level: intermediate 51d4afb720SToby Isaac 5287497f52SBarry Smith Note: 5387497f52SBarry Smith A `PetscDTNodeType` can be paired with a `PetscBool` to indicate whether 5487497f52SBarry Smith the nodes include endpoints or not, and in the case of `PETSCDT_GAUSSJACOBI` 55d4afb720SToby Isaac with exponents for the weight function. 56d4afb720SToby Isaac 57*16a05f60SBarry Smith .seealso: `PetscQuadrature` 58d4afb720SToby Isaac E*/ 599371c9d4SSatish Balay typedef enum { 609371c9d4SSatish Balay PETSCDTNODES_DEFAULT = -1, 619371c9d4SSatish Balay PETSCDTNODES_GAUSSJACOBI, 629371c9d4SSatish Balay PETSCDTNODES_EQUISPACED, 639371c9d4SSatish Balay PETSCDTNODES_TANHSINH 649371c9d4SSatish Balay } PetscDTNodeType; 65d4afb720SToby Isaac 66d3c69ad0SToby Isaac PETSC_EXTERN const char *const *const PetscDTNodeTypes; 67d3c69ad0SToby Isaac 68d3c69ad0SToby Isaac /*E 69d3c69ad0SToby Isaac PetscDTSimplexQuadratureType - A description of classes of quadrature rules for simplices 70d3c69ad0SToby Isaac 71*16a05f60SBarry Smith Values: 72*16a05f60SBarry Smith + `PETSCDTSIMPLEXQUAD_DEFAULT` - Quadrature rule chosen by PETSc 73*16a05f60SBarry Smith . `PETSCDTSIMPLEXQUAD_CONIC` - Quadrature rules constructed as 74d3c69ad0SToby Isaac conically-warped tensor products of 1D 75d3c69ad0SToby Isaac Gauss-Jacobi quadrature rules. These are 76d3c69ad0SToby Isaac explicitly computable in any dimension for any 77d3c69ad0SToby Isaac degree, and the tensor-product structure can be 78d3c69ad0SToby Isaac exploited by sum-factorization methods, but 79d3c69ad0SToby Isaac they are not efficient in terms of nodes per 80d3c69ad0SToby Isaac polynomial degree. 81*16a05f60SBarry Smith - `PETSCDTSIMPLEXQUAD_MINSYM` - Quadrature rules that are fully symmetric 82d3c69ad0SToby Isaac (symmetries of the simplex preserve the nodes 83d3c69ad0SToby Isaac and weights) with minimal (or near minimal) 84d3c69ad0SToby Isaac number of nodes. In dimensions higher than 1 85d3c69ad0SToby Isaac these are not simple to compute, so lookup 86d3c69ad0SToby Isaac tables are used. 87d3c69ad0SToby Isaac 88*16a05f60SBarry Smith Level: intermediate 89*16a05f60SBarry Smith 90*16a05f60SBarry Smith .seealso: `PetscQuadrature`, `PetscDTSimplexQuadrature()` 91d3c69ad0SToby Isaac E*/ 929371c9d4SSatish Balay typedef enum { 939371c9d4SSatish Balay PETSCDTSIMPLEXQUAD_DEFAULT = -1, 949371c9d4SSatish Balay PETSCDTSIMPLEXQUAD_CONIC = 0, 959371c9d4SSatish Balay PETSCDTSIMPLEXQUAD_MINSYM 969371c9d4SSatish Balay } PetscDTSimplexQuadratureType; 97d3c69ad0SToby Isaac 98d3c69ad0SToby Isaac PETSC_EXTERN const char *const *const PetscDTSimplexQuadratureTypes; 99d4afb720SToby Isaac 10021454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *); 101c9638911SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *); 1024366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetCellType(PetscQuadrature, DMPolytopeType *); 1034366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetCellType(PetscQuadrature, DMPolytopeType); 104bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt *); 105bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt); 106a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt *); 107a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt); 1084f9ab2b4SJed Brown PETSC_EXTERN PetscErrorCode PetscQuadratureEqual(PetscQuadrature, PetscQuadrature, PetscBool *); 109a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt *, PetscInt *, PetscInt *, const PetscReal *[], const PetscReal *[]); 110a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal[], const PetscReal[]); 11121454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer); 11221454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *); 113a0845e3aSMatthew G. Knepley 1142df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTensorQuadratureCreate(PetscQuadrature, PetscQuadrature, PetscQuadrature *); 11589710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *); 1164366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureComputePermutations(PetscQuadrature, PetscInt *, IS *[]); 11789710940SMatthew G. Knepley 118907761f8SToby Isaac PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *); 119907761f8SToby Isaac 12037045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *); 121fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiNorm(PetscReal, PetscReal, PetscInt, PetscReal *); 12294e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt, PetscReal, PetscReal, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *); 123fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEvalJet(PetscReal, PetscReal, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]); 124fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPKDEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]); 125d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedSize(PetscInt, PetscInt, PetscInt, PetscInt *); 126d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscInt, PetscReal[]); 12737045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt, PetscReal, PetscReal, PetscReal *, PetscReal *); 12894e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *); 12994e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *); 130916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt, PetscGaussLobattoLegendreCreateType, PetscReal *, PetscReal *); 131194825f6SJed Brown PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *); 132a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *); 133e6a796c3SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *); 134d3c69ad0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTSimplexQuadrature(PetscInt, PetscInt, PetscDTSimplexQuadratureType, PetscQuadrature *); 1354366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTCreateDefaultQuadrature(DMPolytopeType, PetscInt, PetscQuadrature *, PetscQuadrature *); 13637045ce4SJed Brown 137b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *); 138d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *); 139d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *); 140b3c0f97bSTom Klotz 141916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *); 142916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 143916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 144916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 145916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 146916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 147916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 148916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 149916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 150916e780bShannah_mairs 1511a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 1521a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 1531a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 1541a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *); 1551a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *); 1561a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 1571a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *); 158dda711d0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]); 1591a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 1601a989b97SToby Isaac 161d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt, PetscInt, const PetscInt[], PetscInt *); 162d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt, PetscInt, PetscInt, PetscInt[]); 163fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGradedOrderToIndex(PetscInt, const PetscInt[], PetscInt *); 164fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToGradedOrder(PetscInt, PetscInt, PetscInt[]); 165d4afb720SToby Isaac 166fad4db65SToby Isaac #if defined(PETSC_USE_64BIT_INDICES) 167fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 20 168fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX 61 169fad4db65SToby Isaac #else 170fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 12 171fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX 29 172fad4db65SToby Isaac #endif 173fad4db65SToby Isaac 174fad4db65SToby Isaac /*MC 175fad4db65SToby Isaac PetscDTFactorial - Approximate n! as a real number 176fad4db65SToby Isaac 1774165533cSJose E. Roman Input Parameter: 178fad4db65SToby Isaac . n - a non-negative integer 179fad4db65SToby Isaac 1804165533cSJose E. Roman Output Parameter: 181fad4db65SToby Isaac . factorial - n! 182fad4db65SToby Isaac 183fad4db65SToby Isaac Level: beginner 184*16a05f60SBarry Smith 185*16a05f60SBarry Smith .seealso: `PetscDTFactorialInt()`, `PetscDTBinomialInt()`, `PetscDTBinomial()` 186fad4db65SToby Isaac M*/ 187d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial) 188d71ae5a4SJacob Faibussowitsch { 189fad4db65SToby Isaac PetscReal f = 1.0; 190fad4db65SToby Isaac 191fad4db65SToby Isaac PetscFunctionBegin; 192e2ab39ccSLisandro Dalcin *factorial = -1.0; 19363a3b9bcSJacob Faibussowitsch PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %" PetscInt_FMT, n); 1945f80ce2aSJacob Faibussowitsch for (PetscInt i = 1; i < n + 1; ++i) f *= (PetscReal)i; 195fad4db65SToby Isaac *factorial = f; 1963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 197fad4db65SToby Isaac } 198fad4db65SToby Isaac 199fad4db65SToby Isaac /*MC 200fad4db65SToby Isaac PetscDTFactorialInt - Compute n! as an integer 201fad4db65SToby Isaac 2024165533cSJose E. Roman Input Parameter: 203fad4db65SToby Isaac . n - a non-negative integer 204fad4db65SToby Isaac 2054165533cSJose E. Roman Output Parameter: 206fad4db65SToby Isaac . factorial - n! 207fad4db65SToby Isaac 208fad4db65SToby Isaac Level: beginner 209fad4db65SToby Isaac 210*16a05f60SBarry Smith Note: 211*16a05f60SBarry 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. 212*16a05f60SBarry Smith 213*16a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTBinomialInt()`, `PetscDTBinomial()` 214fad4db65SToby Isaac M*/ 215d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial) 216d71ae5a4SJacob Faibussowitsch { 217fad4db65SToby Isaac PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600}; 218fad4db65SToby Isaac 21928222859SToby Isaac PetscFunctionBegin; 22028222859SToby Isaac *factorial = -1; 22163a3b9bcSJacob 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); 222fad4db65SToby Isaac if (n <= 12) { 223fad4db65SToby Isaac *factorial = facLookup[n]; 224fad4db65SToby Isaac } else { 225fad4db65SToby Isaac PetscInt f = facLookup[12]; 226fad4db65SToby Isaac PetscInt i; 227fad4db65SToby Isaac 228fad4db65SToby Isaac for (i = 13; i < n + 1; ++i) f *= i; 229fad4db65SToby Isaac *factorial = f; 230fad4db65SToby Isaac } 2313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 232fad4db65SToby Isaac } 233fad4db65SToby Isaac 234fad4db65SToby Isaac /*MC 235fad4db65SToby Isaac PetscDTBinomial - Approximate the binomial coefficient "n choose k" 236fad4db65SToby Isaac 2374165533cSJose E. Roman Input Parameters: 238fad4db65SToby Isaac + n - a non-negative integer 239fad4db65SToby Isaac - k - an integer between 0 and n, inclusive 240fad4db65SToby Isaac 2414165533cSJose E. Roman Output Parameter: 242fad4db65SToby Isaac . binomial - approximation of the binomial coefficient n choose k 243fad4db65SToby Isaac 244fad4db65SToby Isaac Level: beginner 245*16a05f60SBarry Smith 246*16a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()` 247fad4db65SToby Isaac M*/ 248d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial) 249d71ae5a4SJacob Faibussowitsch { 2501a989b97SToby Isaac PetscFunctionBeginHot; 251e2ab39ccSLisandro Dalcin *binomial = -1.0; 25263a3b9bcSJacob 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); 2531a989b97SToby Isaac if (n <= 3) { 2549371c9d4SSatish Balay PetscInt binomLookup[4][4] = { 2559371c9d4SSatish Balay {1, 0, 0, 0}, 2569371c9d4SSatish Balay {1, 1, 0, 0}, 2579371c9d4SSatish Balay {1, 2, 1, 0}, 2589371c9d4SSatish Balay {1, 3, 3, 1} 2599371c9d4SSatish Balay }; 2601a989b97SToby Isaac 261e2ab39ccSLisandro Dalcin *binomial = (PetscReal)binomLookup[n][k]; 2621a989b97SToby Isaac } else { 263e2ab39ccSLisandro Dalcin PetscReal binom = 1.0; 2641a989b97SToby Isaac 2651a989b97SToby Isaac k = PetscMin(k, n - k); 2665f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1); 2671a989b97SToby Isaac *binomial = binom; 2681a989b97SToby Isaac } 2693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2701a989b97SToby Isaac } 2711a989b97SToby Isaac 272fad4db65SToby Isaac /*MC 273fad4db65SToby Isaac PetscDTBinomialInt - Compute the binomial coefficient "n choose k" 274fad4db65SToby Isaac 27597bb3fdcSJose E. Roman Input Parameters: 276fad4db65SToby Isaac + n - a non-negative integer 277fad4db65SToby Isaac - k - an integer between 0 and n, inclusive 278fad4db65SToby Isaac 27997bb3fdcSJose E. Roman Output Parameter: 280fad4db65SToby Isaac . binomial - the binomial coefficient n choose k 281fad4db65SToby Isaac 282fad4db65SToby Isaac Level: beginner 283*16a05f60SBarry Smith 284*16a05f60SBarry Smith Note: 285*16a05f60SBarry Smith This is limited by integers that can be represented by `PetscInt`. 286*16a05f60SBarry Smith 287*16a05f60SBarry Smith Use `PetscDTBinomial()` for real number approximations of larger values 288*16a05f60SBarry Smith 289*16a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTEnumPerm()` 290fad4db65SToby Isaac M*/ 291d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial) 292d71ae5a4SJacob Faibussowitsch { 29328222859SToby Isaac PetscInt bin; 29428222859SToby Isaac 29528222859SToby Isaac PetscFunctionBegin; 29628222859SToby Isaac *binomial = -1; 29763a3b9bcSJacob 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); 29863a3b9bcSJacob 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); 299fad4db65SToby Isaac if (n <= 3) { 3009371c9d4SSatish Balay PetscInt binomLookup[4][4] = { 3019371c9d4SSatish Balay {1, 0, 0, 0}, 3029371c9d4SSatish Balay {1, 1, 0, 0}, 3039371c9d4SSatish Balay {1, 2, 1, 0}, 3049371c9d4SSatish Balay {1, 3, 3, 1} 3059371c9d4SSatish Balay }; 306fad4db65SToby Isaac 30728222859SToby Isaac bin = binomLookup[n][k]; 308fad4db65SToby Isaac } else { 309fad4db65SToby Isaac PetscInt binom = 1; 310fad4db65SToby Isaac 311fad4db65SToby Isaac k = PetscMin(k, n - k); 3125f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1); 31328222859SToby Isaac bin = binom; 314fad4db65SToby Isaac } 31528222859SToby Isaac *binomial = bin; 3163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 317fad4db65SToby Isaac } 318fad4db65SToby Isaac 319fad4db65SToby Isaac /*MC 320*16a05f60SBarry Smith PetscDTEnumPerm - Get a permutation of `n` integers from its encoding into the integers [0, n!) as a sequence of swaps. 321fad4db65SToby Isaac 3224165533cSJose E. Roman Input Parameters: 323fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 3248cd1e013SToby Isaac - k - an integer in [0, n!) 325fad4db65SToby Isaac 3264165533cSJose E. Roman Output Parameters: 327fad4db65SToby Isaac + perm - the permuted list of the integers [0, ..., n-1] 3282fe279fdSBarry Smith - isOdd - if not `NULL`, returns whether the permutation used an even or odd number of swaps. 329fad4db65SToby Isaac 330*16a05f60SBarry Smith Level: intermediate 331fad4db65SToby Isaac 332*16a05f60SBarry Smith Notes: 333*16a05f60SBarry Smith A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation, 334*16a05f60SBarry Smith by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in 335*16a05f60SBarry Smith some position j >= i. This swap is encoded as the difference (j - i). The difference d_i at step i is less than 336*16a05f60SBarry Smith (n - i). This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number 337*16a05f60SBarry Smith (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}. 338*16a05f60SBarry Smith 339*16a05f60SBarry 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. 340*16a05f60SBarry Smith 341*16a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTPermIndex()` 342fad4db65SToby Isaac M*/ 343d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd) 344d71ae5a4SJacob Faibussowitsch { 3451a989b97SToby Isaac PetscInt odd = 0; 3461a989b97SToby Isaac PetscInt i; 347fad4db65SToby Isaac PetscInt work[PETSC_FACTORIAL_MAX]; 348fad4db65SToby Isaac PetscInt *w; 3491a989b97SToby Isaac 35028222859SToby Isaac PetscFunctionBegin; 35128222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 35263a3b9bcSJacob 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); 353fad4db65SToby Isaac w = &work[n - 2]; 3541a989b97SToby Isaac for (i = 2; i <= n; i++) { 3551a989b97SToby Isaac *(w--) = k % i; 3561a989b97SToby Isaac k /= i; 3571a989b97SToby Isaac } 3581a989b97SToby Isaac for (i = 0; i < n; i++) perm[i] = i; 3591a989b97SToby Isaac for (i = 0; i < n - 1; i++) { 3601a989b97SToby Isaac PetscInt s = work[i]; 3611a989b97SToby Isaac PetscInt swap = perm[i]; 3621a989b97SToby Isaac 3631a989b97SToby Isaac perm[i] = perm[i + s]; 3641a989b97SToby Isaac perm[i + s] = swap; 3651a989b97SToby Isaac odd ^= (!!s); 3661a989b97SToby Isaac } 3671a989b97SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 3683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3691a989b97SToby Isaac } 3701a989b97SToby Isaac 371fad4db65SToby Isaac /*MC 372*16a05f60SBarry Smith PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!). This inverts `PetscDTEnumPerm()`. 3738cd1e013SToby Isaac 3744165533cSJose E. Roman Input Parameters: 3758cd1e013SToby Isaac + n - a non-negative integer (see note about limits below) 3768cd1e013SToby Isaac - perm - the permuted list of the integers [0, ..., n-1] 3778cd1e013SToby Isaac 3784165533cSJose E. Roman Output Parameters: 3798cd1e013SToby Isaac + k - an integer in [0, n!) 3802fe279fdSBarry Smith - isOdd - if not `NULL`, returns whether the permutation used an even or odd number of swaps. 3818cd1e013SToby Isaac 3828cd1e013SToby Isaac Level: beginner 383*16a05f60SBarry Smith 384*16a05f60SBarry Smith Note: 385*16a05f60SBarry 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. 386*16a05f60SBarry Smith 387*16a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()` 3888cd1e013SToby Isaac M*/ 389d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd) 390d71ae5a4SJacob Faibussowitsch { 3918cd1e013SToby Isaac PetscInt odd = 0; 3928cd1e013SToby Isaac PetscInt i, idx; 3938cd1e013SToby Isaac PetscInt work[PETSC_FACTORIAL_MAX]; 3948cd1e013SToby Isaac PetscInt iwork[PETSC_FACTORIAL_MAX]; 3958cd1e013SToby Isaac 3968cd1e013SToby Isaac PetscFunctionBeginHot; 39728222859SToby Isaac *k = -1; 39828222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 39963a3b9bcSJacob 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); 4008cd1e013SToby Isaac for (i = 0; i < n; i++) work[i] = i; /* partial permutation */ 4018cd1e013SToby Isaac for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */ 4028cd1e013SToby Isaac for (idx = 0, i = 0; i < n - 1; i++) { 4038cd1e013SToby Isaac PetscInt j = perm[i]; 4048cd1e013SToby Isaac PetscInt icur = work[i]; 4058cd1e013SToby Isaac PetscInt jloc = iwork[j]; 4068cd1e013SToby Isaac PetscInt diff = jloc - i; 4078cd1e013SToby Isaac 4088cd1e013SToby Isaac idx = idx * (n - i) + diff; 4098cd1e013SToby Isaac /* swap (i, jloc) */ 4108cd1e013SToby Isaac work[i] = j; 4118cd1e013SToby Isaac work[jloc] = icur; 4128cd1e013SToby Isaac iwork[j] = i; 4138cd1e013SToby Isaac iwork[icur] = jloc; 4148cd1e013SToby Isaac odd ^= (!!diff); 4158cd1e013SToby Isaac } 4168cd1e013SToby Isaac *k = idx; 4178cd1e013SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 4183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4198cd1e013SToby Isaac } 4208cd1e013SToby Isaac 4218cd1e013SToby Isaac /*MC 422fad4db65SToby Isaac PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k). 423fad4db65SToby Isaac The encoding is in lexicographic order. 424fad4db65SToby Isaac 4254165533cSJose E. Roman Input Parameters: 426fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 427fad4db65SToby Isaac . k - an integer in [0, n] 428fad4db65SToby Isaac - j - an index in [0, n choose k) 429fad4db65SToby Isaac 4304165533cSJose E. Roman Output Parameter: 431fad4db65SToby Isaac . subset - the jth subset of size k of the integers [0, ..., n - 1] 432fad4db65SToby Isaac 433fad4db65SToby Isaac Level: beginner 434fad4db65SToby Isaac 435*16a05f60SBarry Smith Note: 436*16a05f60SBarry Smith Limited by arguments such that `n` choose `k` can be represented by `PetscInt` 437*16a05f60SBarry Smith 438*16a05f60SBarry Smith .seealso: `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()` 439fad4db65SToby Isaac M*/ 440d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset) 441d71ae5a4SJacob Faibussowitsch { 4425f80ce2aSJacob Faibussowitsch PetscInt Nk; 4431a989b97SToby Isaac 4441a989b97SToby Isaac PetscFunctionBeginHot; 4459566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 4465f80ce2aSJacob Faibussowitsch for (PetscInt i = 0, l = 0; i < n && l < k; i++) { 4471a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 4481a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 4491a989b97SToby Isaac 4501a989b97SToby Isaac if (j < Nminuskminus) { 4511a989b97SToby Isaac subset[l++] = i; 4521a989b97SToby Isaac Nk = Nminuskminus; 4531a989b97SToby Isaac } else { 4541a989b97SToby Isaac j -= Nminuskminus; 4551a989b97SToby Isaac Nk = Nminusk; 4561a989b97SToby Isaac } 4571a989b97SToby Isaac } 4583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4591a989b97SToby Isaac } 4601a989b97SToby Isaac 461fad4db65SToby Isaac /*MC 46287497f52SBarry 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. 46387497f52SBarry Smith This is the inverse of `PetscDTEnumSubset`. 464fad4db65SToby Isaac 4654165533cSJose E. Roman Input Parameters: 466fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 467fad4db65SToby Isaac . k - an integer in [0, n] 468fad4db65SToby Isaac - subset - an ordered subset of the integers [0, ..., n - 1] 469fad4db65SToby Isaac 4704165533cSJose E. Roman Output Parameter: 471fad4db65SToby Isaac . index - the rank of the subset in lexicographic order 472fad4db65SToby Isaac 473fad4db65SToby Isaac Level: beginner 474fad4db65SToby Isaac 475*16a05f60SBarry Smith Note: 476*16a05f60SBarry Smith Limited by arguments such that `n` choose `k` can be represented by `PetscInt` 477*16a05f60SBarry Smith 478*16a05f60SBarry Smith .seealso: `PetscDTEnumSubset()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()` 479fad4db65SToby Isaac M*/ 480d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index) 481d71ae5a4SJacob Faibussowitsch { 4825f80ce2aSJacob Faibussowitsch PetscInt j = 0, Nk; 4831a989b97SToby Isaac 48428222859SToby Isaac PetscFunctionBegin; 48528222859SToby Isaac *index = -1; 4869566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 4875f80ce2aSJacob Faibussowitsch for (PetscInt i = 0, l = 0; i < n && l < k; i++) { 4881a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 4891a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 4901a989b97SToby Isaac 4911a989b97SToby Isaac if (subset[l] == i) { 4921a989b97SToby Isaac l++; 4931a989b97SToby Isaac Nk = Nminuskminus; 4941a989b97SToby Isaac } else { 4951a989b97SToby Isaac j += Nminuskminus; 4961a989b97SToby Isaac Nk = Nminusk; 4971a989b97SToby Isaac } 4981a989b97SToby Isaac } 4991a989b97SToby Isaac *index = j; 5003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5011a989b97SToby Isaac } 5021a989b97SToby Isaac 503fad4db65SToby Isaac /*MC 50428222859SToby 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. 505fad4db65SToby Isaac 5064165533cSJose E. Roman Input Parameters: 507fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 508fad4db65SToby Isaac . k - an integer in [0, n] 509fad4db65SToby Isaac - j - an index in [0, n choose k) 510fad4db65SToby Isaac 5114165533cSJose E. Roman Output Parameters: 512fad4db65SToby Isaac + perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set. 5132fe279fdSBarry Smith - isOdd - if not `NULL`, return whether perm is an even or odd permutation. 514fad4db65SToby Isaac 515fad4db65SToby Isaac Level: beginner 516fad4db65SToby Isaac 517*16a05f60SBarry Smith Note: 518*16a05f60SBarry Smith Limited by arguments such that `n` choose `k` can be represented by `PetscInt` 519*16a05f60SBarry Smith 520*16a05f60SBarry Smith .seealso: `PetscDTEnumSubset()`, `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, 521*16a05f60SBarry Smith `PetscDTPermIndex()` 522fad4db65SToby Isaac M*/ 523d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd) 524d71ae5a4SJacob Faibussowitsch { 5255f80ce2aSJacob Faibussowitsch PetscInt i, l, m, Nk, odd = 0; 5265f80ce2aSJacob Faibussowitsch PetscInt *subcomp = perm + k; 5271a989b97SToby Isaac 52828222859SToby Isaac PetscFunctionBegin; 52928222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 5309566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 5311a989b97SToby Isaac for (i = 0, l = 0, m = 0; i < n && l < k; i++) { 5321a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 5331a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 5341a989b97SToby Isaac 5351a989b97SToby Isaac if (j < Nminuskminus) { 536fad4db65SToby Isaac perm[l++] = i; 5371a989b97SToby Isaac Nk = Nminuskminus; 5381a989b97SToby Isaac } else { 5391a989b97SToby Isaac subcomp[m++] = i; 5401a989b97SToby Isaac j -= Nminuskminus; 5411a989b97SToby Isaac odd ^= ((k - l) & 1); 5421a989b97SToby Isaac Nk = Nminusk; 5431a989b97SToby Isaac } 5441a989b97SToby Isaac } 5455f80ce2aSJacob Faibussowitsch for (; i < n; i++) subcomp[m++] = i; 5461a989b97SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 5473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5481a989b97SToby Isaac } 5491a989b97SToby Isaac 550ef0bb6c7SMatthew G. Knepley struct _p_PetscTabulation { 551a5b23f4aSJose E. Roman PetscInt K; /* Indicates a k-jet, namely tabulated derivatives up to order k */ 55219815104SMartin Diehl PetscInt Nr; /* The number of tabulation replicas (often 1) */ 553ef0bb6c7SMatthew G. Knepley PetscInt Np; /* The number of tabulation points in a replica */ 554ef0bb6c7SMatthew G. Knepley PetscInt Nb; /* The number of functions tabulated */ 555ef0bb6c7SMatthew G. Knepley PetscInt Nc; /* The number of function components */ 556ef0bb6c7SMatthew G. Knepley PetscInt cdim; /* The coordinate dimension */ 557ef0bb6c7SMatthew G. Knepley PetscReal **T; /* The tabulation T[K] of functions and their derivatives 558ef0bb6c7SMatthew G. Knepley T[0] = B[Nr*Np][Nb][Nc]: The basis function values at quadrature points 559ef0bb6c7SMatthew G. Knepley T[1] = D[Nr*Np][Nb][Nc][cdim]: The basis function derivatives at quadrature points 560ef0bb6c7SMatthew G. Knepley T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */ 561ef0bb6c7SMatthew G. Knepley }; 562ef0bb6c7SMatthew G. Knepley typedef struct _p_PetscTabulation *PetscTabulation; 563ef0bb6c7SMatthew G. Knepley 564d6685f55SMatthew G. Knepley typedef PetscErrorCode (*PetscProbFunc)(const PetscReal[], const PetscReal[], PetscReal[]); 565d6685f55SMatthew G. Knepley 5669371c9d4SSatish Balay typedef enum { 5679371c9d4SSatish Balay DTPROB_DENSITY_CONSTANT, 5689371c9d4SSatish Balay DTPROB_DENSITY_GAUSSIAN, 5699371c9d4SSatish Balay DTPROB_DENSITY_MAXWELL_BOLTZMANN, 5709371c9d4SSatish Balay DTPROB_NUM_DENSITY 5719371c9d4SSatish Balay } DTProbDensityType; 572d6685f55SMatthew G. Knepley PETSC_EXTERN const char *const DTProbDensityTypes[]; 573d6685f55SMatthew G. Knepley 574d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]); 575d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]); 576d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]); 577d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]); 578d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]); 579d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]); 580d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 581d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 582d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 583d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]); 584d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]); 585ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]); 586ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]); 587d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 588d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 589d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 590ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 591ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 592ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 593ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 594ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 595ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 596d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFunc *, PetscProbFunc *, PetscProbFunc *); 597d6685f55SMatthew G. Knepley 598d6685f55SMatthew G. Knepley #include <petscvec.h> 599d6685f55SMatthew G. Knepley 600d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFunc, PetscReal *); 601d6685f55SMatthew G. Knepley 60237045ce4SJed Brown #endif 603