137045ce4SJed Brown /* 237045ce4SJed Brown Common tools for constructing discretizations 337045ce4SJed Brown */ 426bd1501SBarry Smith #if !defined(PETSCDT_H) 526bd1501SBarry Smith #define PETSCDT_H 637045ce4SJed Brown 737045ce4SJed Brown #include <petscsys.h> 837045ce4SJed Brown 9ac09b921SBarry Smith /* SUBMANSEC = DT */ 10ac09b921SBarry Smith 112cd22861SMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCQUADRATURE_CLASSID; 122cd22861SMatthew G. Knepley 1321454ff5SMatthew G. Knepley /*S 1421454ff5SMatthew G. Knepley PetscQuadrature - Quadrature rule for integration. 1521454ff5SMatthew G. Knepley 16329bbf4eSMatthew G. Knepley Level: beginner 1721454ff5SMatthew G. Knepley 18db781477SPatrick Sanan .seealso: `PetscQuadratureCreate()`, `PetscQuadratureDestroy()` 1921454ff5SMatthew G. Knepley S*/ 2021454ff5SMatthew G. Knepley typedef struct _p_PetscQuadrature *PetscQuadrature; 2121454ff5SMatthew G. Knepley 228272889dSSatish Balay /*E 23916e780bShannah_mairs PetscGaussLobattoLegendreCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights 248272889dSSatish Balay 258272889dSSatish Balay Level: intermediate 268272889dSSatish Balay 27*87497f52SBarry Smith $ `PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA` - compute the nodes via linear algebra 28*87497f52SBarry Smith $ `PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON` - compute the nodes by solving a nonlinear equation with Newton's method 298272889dSSatish Balay 308272889dSSatish Balay E*/ 31f2e8fe4dShannah_mairs typedef enum {PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON} PetscGaussLobattoLegendreCreateType; 328272889dSSatish Balay 33d4afb720SToby Isaac /*E 34d4afb720SToby Isaac PetscDTNodeType - A description of strategies for generating nodes (both 35d4afb720SToby Isaac quadrature nodes and nodes for Lagrange polynomials) 36d4afb720SToby Isaac 37d4afb720SToby Isaac Level: intermediate 38d4afb720SToby Isaac 39*87497f52SBarry Smith $ `PETSCDTNODES_DEFAULT` - Nodes chosen by PETSc 40*87497f52SBarry Smith $ `PETSCDTNODES_GAUSSJACOBI` - Nodes at either Gauss-Jacobi or Gauss-Lobatto-Jacobi quadrature points 41*87497f52SBarry Smith $ `PETSCDTNODES_EQUISPACED` - Nodes equispaced either including the endpoints or excluding them 42*87497f52SBarry Smith $ `PETSCDTNODES_TANHSINH` - Nodes at Tanh-Sinh quadrature points 43d4afb720SToby Isaac 44*87497f52SBarry Smith Note: 45*87497f52SBarry Smith A `PetscDTNodeType` can be paired with a `PetscBool` to indicate whether 46*87497f52SBarry Smith the nodes include endpoints or not, and in the case of `PETSCDT_GAUSSJACOBI` 47d4afb720SToby Isaac with exponents for the weight function. 48d4afb720SToby Isaac 49d4afb720SToby Isaac E*/ 50d4afb720SToby Isaac typedef enum {PETSCDTNODES_DEFAULT=-1, PETSCDTNODES_GAUSSJACOBI, PETSCDTNODES_EQUISPACED, PETSCDTNODES_TANHSINH} PetscDTNodeType; 51d4afb720SToby Isaac 52d3c69ad0SToby Isaac PETSC_EXTERN const char *const*const PetscDTNodeTypes; 53d3c69ad0SToby Isaac 54d3c69ad0SToby Isaac /*E 55d3c69ad0SToby Isaac PetscDTSimplexQuadratureType - A description of classes of quadrature rules for simplices 56d3c69ad0SToby Isaac 57d3c69ad0SToby Isaac Level: intermediate 58d3c69ad0SToby Isaac 59*87497f52SBarry Smith $ `PETSCDTSIMPLEXQUAD_DEFAULT` - Quadrature rule chosen by PETSc 60*87497f52SBarry Smith $ `PETSCDTSIMPLEXQUAD_CONIC` - Quadrature rules constructed as 61d3c69ad0SToby Isaac conically-warped tensor products of 1D 62d3c69ad0SToby Isaac Gauss-Jacobi quadrature rules. These are 63d3c69ad0SToby Isaac explicitly computable in any dimension for any 64d3c69ad0SToby Isaac degree, and the tensor-product structure can be 65d3c69ad0SToby Isaac exploited by sum-factorization methods, but 66d3c69ad0SToby Isaac they are not efficient in terms of nodes per 67d3c69ad0SToby Isaac polynomial degree. 68*87497f52SBarry Smith $ `PETSCDTSIMPLEXQUAD_MINSYM` - Quadrature rules that are fully symmetric 69d3c69ad0SToby Isaac (symmetries of the simplex preserve the nodes 70d3c69ad0SToby Isaac and weights) with minimal (or near minimal) 71d3c69ad0SToby Isaac number of nodes. In dimensions higher than 1 72d3c69ad0SToby Isaac these are not simple to compute, so lookup 73d3c69ad0SToby Isaac tables are used. 74d3c69ad0SToby Isaac 75d3c69ad0SToby Isaac .seealso: `PetscDTSimplexQuadrature()` 76d3c69ad0SToby Isaac E*/ 77d3c69ad0SToby Isaac typedef enum {PETSCDTSIMPLEXQUAD_DEFAULT=-1, PETSCDTSIMPLEXQUAD_CONIC=0, PETSCDTSIMPLEXQUAD_MINSYM} PetscDTSimplexQuadratureType; 78d3c69ad0SToby Isaac 79d3c69ad0SToby Isaac PETSC_EXTERN const char *const*const PetscDTSimplexQuadratureTypes; 80d4afb720SToby Isaac 8121454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *); 82c9638911SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *); 83bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt*); 84bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt); 85a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt*); 86a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt); 874f9ab2b4SJed Brown PETSC_EXTERN PetscErrorCode PetscQuadratureEqual(PetscQuadrature, PetscQuadrature, PetscBool*); 88a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt*, PetscInt*, PetscInt*, const PetscReal *[], const PetscReal *[]); 89a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal [], const PetscReal []); 9021454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer); 9121454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *); 92a0845e3aSMatthew G. Knepley 932df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTensorQuadratureCreate(PetscQuadrature, PetscQuadrature, PetscQuadrature *); 9489710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *); 9589710940SMatthew G. Knepley 96907761f8SToby Isaac PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *); 97907761f8SToby Isaac 9837045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*); 99fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiNorm(PetscReal,PetscReal,PetscInt,PetscReal *); 10094e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt,PetscReal,PetscReal,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*); 101fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEvalJet(PetscReal,PetscReal,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscReal[]); 102fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPKDEvalJet(PetscInt,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscReal[]); 103d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedSize(PetscInt,PetscInt,PetscInt,PetscInt*); 104d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscInt,PetscReal[]); 10537045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt,PetscReal,PetscReal,PetscReal*,PetscReal*); 10694e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*,PetscReal*); 10794e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*,PetscReal*); 108916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt,PetscGaussLobattoLegendreCreateType,PetscReal*,PetscReal*); 109194825f6SJed Brown PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt,PetscInt,const PetscReal*,PetscInt,const PetscReal*,PetscReal*); 110a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*); 111e6a796c3SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*); 112d3c69ad0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTSimplexQuadrature(PetscInt,PetscInt,PetscDTSimplexQuadratureType,PetscQuadrature*); 11337045ce4SJed Brown 114b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *); 115d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *); 116d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *); 117b3c0f97bSTom Klotz 118916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *); 119916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 120916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 121916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 122916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 123916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 124916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 125916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 126916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 127916e780bShannah_mairs 1281a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 1291a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 1301a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 1311a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *); 1321a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *); 1331a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 1341a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *); 135dda711d0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]); 1361a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 1371a989b97SToby Isaac 138d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt,PetscInt,const PetscInt[],PetscInt*); 139d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt,PetscInt,PetscInt,PetscInt[]); 140fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGradedOrderToIndex(PetscInt,const PetscInt[],PetscInt*); 141fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToGradedOrder(PetscInt,PetscInt,PetscInt[]); 142d4afb720SToby Isaac 143fad4db65SToby Isaac #if defined(PETSC_USE_64BIT_INDICES) 144fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 20 145fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX 61 146fad4db65SToby Isaac #else 147fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 12 148fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX 29 149fad4db65SToby Isaac #endif 150fad4db65SToby Isaac 151fad4db65SToby Isaac /*MC 152fad4db65SToby Isaac PetscDTFactorial - Approximate n! as a real number 153fad4db65SToby Isaac 1544165533cSJose E. Roman Input Parameter: 155fad4db65SToby Isaac . n - a non-negative integer 156fad4db65SToby Isaac 1574165533cSJose E. Roman Output Parameter: 158fad4db65SToby Isaac . factorial - n! 159fad4db65SToby Isaac 160fad4db65SToby Isaac Level: beginner 161fad4db65SToby Isaac M*/ 1629fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial) 163fad4db65SToby Isaac { 164fad4db65SToby Isaac PetscReal f = 1.0; 165fad4db65SToby Isaac 166fad4db65SToby Isaac PetscFunctionBegin; 167e2ab39ccSLisandro Dalcin *factorial = -1.0; 16863a3b9bcSJacob Faibussowitsch PetscCheck(n >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %" PetscInt_FMT, n); 1695f80ce2aSJacob Faibussowitsch for (PetscInt i = 1; i < n+1; ++i) f *= (PetscReal)i; 170fad4db65SToby Isaac *factorial = f; 171fad4db65SToby Isaac PetscFunctionReturn(0); 172fad4db65SToby Isaac } 173fad4db65SToby Isaac 174fad4db65SToby Isaac /*MC 175fad4db65SToby Isaac PetscDTFactorialInt - Compute n! as an integer 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 184fad4db65SToby Isaac 185fad4db65SToby 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. 186fad4db65SToby Isaac M*/ 1879fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial) 188fad4db65SToby Isaac { 189fad4db65SToby Isaac PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600}; 190fad4db65SToby Isaac 19128222859SToby Isaac PetscFunctionBegin; 19228222859SToby Isaac *factorial = -1; 19363a3b9bcSJacob 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); 194fad4db65SToby Isaac if (n <= 12) { 195fad4db65SToby Isaac *factorial = facLookup[n]; 196fad4db65SToby Isaac } else { 197fad4db65SToby Isaac PetscInt f = facLookup[12]; 198fad4db65SToby Isaac PetscInt i; 199fad4db65SToby Isaac 200fad4db65SToby Isaac for (i = 13; i < n+1; ++i) f *= i; 201fad4db65SToby Isaac *factorial = f; 202fad4db65SToby Isaac } 203fad4db65SToby Isaac PetscFunctionReturn(0); 204fad4db65SToby Isaac } 205fad4db65SToby Isaac 206fad4db65SToby Isaac /*MC 207fad4db65SToby Isaac PetscDTBinomial - Approximate the binomial coefficient "n choose k" 208fad4db65SToby Isaac 2094165533cSJose E. Roman Input Parameters: 210fad4db65SToby Isaac + n - a non-negative integer 211fad4db65SToby Isaac - k - an integer between 0 and n, inclusive 212fad4db65SToby Isaac 2134165533cSJose E. Roman Output Parameter: 214fad4db65SToby Isaac . binomial - approximation of the binomial coefficient n choose k 215fad4db65SToby Isaac 216fad4db65SToby Isaac Level: beginner 217fad4db65SToby Isaac M*/ 2189fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial) 2191a989b97SToby Isaac { 2201a989b97SToby Isaac PetscFunctionBeginHot; 221e2ab39ccSLisandro Dalcin *binomial = -1.0; 22263a3b9bcSJacob 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); 2231a989b97SToby Isaac if (n <= 3) { 2241a989b97SToby Isaac PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}}; 2251a989b97SToby Isaac 226e2ab39ccSLisandro Dalcin *binomial = (PetscReal)binomLookup[n][k]; 2271a989b97SToby Isaac } else { 228e2ab39ccSLisandro Dalcin PetscReal binom = 1.0; 2291a989b97SToby Isaac 2301a989b97SToby Isaac k = PetscMin(k, n - k); 2315f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1); 2321a989b97SToby Isaac *binomial = binom; 2331a989b97SToby Isaac } 2341a989b97SToby Isaac PetscFunctionReturn(0); 2351a989b97SToby Isaac } 2361a989b97SToby Isaac 237fad4db65SToby Isaac /*MC 238fad4db65SToby Isaac PetscDTBinomialInt - Compute the binomial coefficient "n choose k" 239fad4db65SToby Isaac 24097bb3fdcSJose E. Roman Input Parameters: 241fad4db65SToby Isaac + n - a non-negative integer 242fad4db65SToby Isaac - k - an integer between 0 and n, inclusive 243fad4db65SToby Isaac 24497bb3fdcSJose E. Roman Output Parameter: 245fad4db65SToby Isaac . binomial - the binomial coefficient n choose k 246fad4db65SToby Isaac 247*87497f52SBarry Smith Note: this is limited by integers that can be represented by `PetscInt` 248fad4db65SToby Isaac 249fad4db65SToby Isaac Level: beginner 250fad4db65SToby Isaac M*/ 2519fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial) 252fad4db65SToby Isaac { 25328222859SToby Isaac PetscInt bin; 25428222859SToby Isaac 25528222859SToby Isaac PetscFunctionBegin; 25628222859SToby Isaac *binomial = -1; 25763a3b9bcSJacob 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); 25863a3b9bcSJacob 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); 259fad4db65SToby Isaac if (n <= 3) { 260fad4db65SToby Isaac PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}}; 261fad4db65SToby Isaac 26228222859SToby Isaac bin = binomLookup[n][k]; 263fad4db65SToby Isaac } else { 264fad4db65SToby Isaac PetscInt binom = 1; 265fad4db65SToby Isaac 266fad4db65SToby Isaac k = PetscMin(k, n - k); 2675f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1); 26828222859SToby Isaac bin = binom; 269fad4db65SToby Isaac } 27028222859SToby Isaac *binomial = bin; 271fad4db65SToby Isaac PetscFunctionReturn(0); 272fad4db65SToby Isaac } 273fad4db65SToby Isaac 274fad4db65SToby Isaac /*MC 275fad4db65SToby Isaac PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps. 276fad4db65SToby Isaac 277fad4db65SToby Isaac A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation, 278fad4db65SToby Isaac by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in 27928222859SToby Isaac some position j >= i. This swap is encoded as the difference (j - i). The difference d_i at step i is less than 28028222859SToby Isaac (n - i). This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number 281fad4db65SToby Isaac (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}. 282fad4db65SToby Isaac 2834165533cSJose E. Roman Input Parameters: 284fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 2858cd1e013SToby Isaac - k - an integer in [0, n!) 286fad4db65SToby Isaac 2874165533cSJose E. Roman Output Parameters: 288fad4db65SToby Isaac + perm - the permuted list of the integers [0, ..., n-1] 2898cd1e013SToby Isaac - isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps. 290fad4db65SToby Isaac 291*87497f52SBarry Smith Note: 292*87497f52SBarry 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. 293fad4db65SToby Isaac 294fad4db65SToby Isaac Level: beginner 295fad4db65SToby Isaac M*/ 2969fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd) 2971a989b97SToby Isaac { 2981a989b97SToby Isaac PetscInt odd = 0; 2991a989b97SToby Isaac PetscInt i; 300fad4db65SToby Isaac PetscInt work[PETSC_FACTORIAL_MAX]; 301fad4db65SToby Isaac PetscInt *w; 3021a989b97SToby Isaac 30328222859SToby Isaac PetscFunctionBegin; 30428222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 30563a3b9bcSJacob 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); 306fad4db65SToby Isaac w = &work[n - 2]; 3071a989b97SToby Isaac for (i = 2; i <= n; i++) { 3081a989b97SToby Isaac *(w--) = k % i; 3091a989b97SToby Isaac k /= i; 3101a989b97SToby Isaac } 3111a989b97SToby Isaac for (i = 0; i < n; i++) perm[i] = i; 3121a989b97SToby Isaac for (i = 0; i < n - 1; i++) { 3131a989b97SToby Isaac PetscInt s = work[i]; 3141a989b97SToby Isaac PetscInt swap = perm[i]; 3151a989b97SToby Isaac 3161a989b97SToby Isaac perm[i] = perm[i + s]; 3171a989b97SToby Isaac perm[i + s] = swap; 3181a989b97SToby Isaac odd ^= (!!s); 3191a989b97SToby Isaac } 3201a989b97SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 3211a989b97SToby Isaac PetscFunctionReturn(0); 3221a989b97SToby Isaac } 3231a989b97SToby Isaac 324fad4db65SToby Isaac /*MC 325*87497f52SBarry Smith PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!). This inverts `PetscDTEnumPerm`. 3268cd1e013SToby Isaac 3274165533cSJose E. Roman Input Parameters: 3288cd1e013SToby Isaac + n - a non-negative integer (see note about limits below) 3298cd1e013SToby Isaac - perm - the permuted list of the integers [0, ..., n-1] 3308cd1e013SToby Isaac 3314165533cSJose E. Roman Output Parameters: 3328cd1e013SToby Isaac + k - an integer in [0, n!) 333f0fc11ceSJed Brown - isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps. 3348cd1e013SToby Isaac 335*87497f52SBarry Smith Note: 336*87497f52SBarry 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. 3378cd1e013SToby Isaac 3388cd1e013SToby Isaac Level: beginner 3398cd1e013SToby Isaac M*/ 3409fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd) 3418cd1e013SToby Isaac { 3428cd1e013SToby Isaac PetscInt odd = 0; 3438cd1e013SToby Isaac PetscInt i, idx; 3448cd1e013SToby Isaac PetscInt work[PETSC_FACTORIAL_MAX]; 3458cd1e013SToby Isaac PetscInt iwork[PETSC_FACTORIAL_MAX]; 3468cd1e013SToby Isaac 3478cd1e013SToby Isaac PetscFunctionBeginHot; 34828222859SToby Isaac *k = -1; 34928222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 35063a3b9bcSJacob 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); 3518cd1e013SToby Isaac for (i = 0; i < n; i++) work[i] = i; /* partial permutation */ 3528cd1e013SToby Isaac for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */ 3538cd1e013SToby Isaac for (idx = 0, i = 0; i < n - 1; i++) { 3548cd1e013SToby Isaac PetscInt j = perm[i]; 3558cd1e013SToby Isaac PetscInt icur = work[i]; 3568cd1e013SToby Isaac PetscInt jloc = iwork[j]; 3578cd1e013SToby Isaac PetscInt diff = jloc - i; 3588cd1e013SToby Isaac 3598cd1e013SToby Isaac idx = idx * (n - i) + diff; 3608cd1e013SToby Isaac /* swap (i, jloc) */ 3618cd1e013SToby Isaac work[i] = j; 3628cd1e013SToby Isaac work[jloc] = icur; 3638cd1e013SToby Isaac iwork[j] = i; 3648cd1e013SToby Isaac iwork[icur] = jloc; 3658cd1e013SToby Isaac odd ^= (!!diff); 3668cd1e013SToby Isaac } 3678cd1e013SToby Isaac *k = idx; 3688cd1e013SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 3698cd1e013SToby Isaac PetscFunctionReturn(0); 3708cd1e013SToby Isaac } 3718cd1e013SToby Isaac 3728cd1e013SToby Isaac /*MC 373fad4db65SToby Isaac PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k). 374fad4db65SToby Isaac The encoding is in lexicographic order. 375fad4db65SToby Isaac 3764165533cSJose E. Roman Input Parameters: 377fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 378fad4db65SToby Isaac . k - an integer in [0, n] 379fad4db65SToby Isaac - j - an index in [0, n choose k) 380fad4db65SToby Isaac 3814165533cSJose E. Roman Output Parameter: 382fad4db65SToby Isaac . subset - the jth subset of size k of the integers [0, ..., n - 1] 383fad4db65SToby Isaac 384*87497f52SBarry Smith Note: 385*87497f52SBarry Smith Limited by arguments such that n choose k can be represented by `PetscInt` 386fad4db65SToby Isaac 387fad4db65SToby Isaac Level: beginner 388fad4db65SToby Isaac 389db781477SPatrick Sanan .seealso: `PetscDTSubsetIndex()` 390fad4db65SToby Isaac M*/ 3919fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset) 3921a989b97SToby Isaac { 3935f80ce2aSJacob Faibussowitsch PetscInt Nk; 3941a989b97SToby Isaac 3951a989b97SToby Isaac PetscFunctionBeginHot; 3969566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 3975f80ce2aSJacob Faibussowitsch for (PetscInt i = 0, l = 0; i < n && l < k; i++) { 3981a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 3991a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 4001a989b97SToby Isaac 4011a989b97SToby Isaac if (j < Nminuskminus) { 4021a989b97SToby Isaac subset[l++] = i; 4031a989b97SToby Isaac Nk = Nminuskminus; 4041a989b97SToby Isaac } else { 4051a989b97SToby Isaac j -= Nminuskminus; 4061a989b97SToby Isaac Nk = Nminusk; 4071a989b97SToby Isaac } 4081a989b97SToby Isaac } 4091a989b97SToby Isaac PetscFunctionReturn(0); 4101a989b97SToby Isaac } 4111a989b97SToby Isaac 412fad4db65SToby Isaac /*MC 413*87497f52SBarry 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. 414*87497f52SBarry Smith This is the inverse of `PetscDTEnumSubset`. 415fad4db65SToby Isaac 4164165533cSJose E. Roman Input Parameters: 417fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 418fad4db65SToby Isaac . k - an integer in [0, n] 419fad4db65SToby Isaac - subset - an ordered subset of the integers [0, ..., n - 1] 420fad4db65SToby Isaac 4214165533cSJose E. Roman Output Parameter: 422fad4db65SToby Isaac . index - the rank of the subset in lexicographic order 423fad4db65SToby Isaac 424*87497f52SBarry Smith Note: 425*87497f52SBarry Smith Limited by arguments such that n choose k can be represented by `PetscInt` 426fad4db65SToby Isaac 427fad4db65SToby Isaac Level: beginner 428fad4db65SToby Isaac 429db781477SPatrick Sanan .seealso: `PetscDTEnumSubset()` 430fad4db65SToby Isaac M*/ 4319fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index) 4321a989b97SToby Isaac { 4335f80ce2aSJacob Faibussowitsch PetscInt j = 0, Nk; 4341a989b97SToby Isaac 43528222859SToby Isaac PetscFunctionBegin; 43628222859SToby Isaac *index = -1; 4379566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 4385f80ce2aSJacob Faibussowitsch for (PetscInt i = 0, l = 0; i < n && l < k; i++) { 4391a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 4401a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 4411a989b97SToby Isaac 4421a989b97SToby Isaac if (subset[l] == i) { 4431a989b97SToby Isaac l++; 4441a989b97SToby Isaac Nk = Nminuskminus; 4451a989b97SToby Isaac } else { 4461a989b97SToby Isaac j += Nminuskminus; 4471a989b97SToby Isaac Nk = Nminusk; 4481a989b97SToby Isaac } 4491a989b97SToby Isaac } 4501a989b97SToby Isaac *index = j; 4511a989b97SToby Isaac PetscFunctionReturn(0); 4521a989b97SToby Isaac } 4531a989b97SToby Isaac 454fad4db65SToby Isaac /*MC 45528222859SToby 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. 456fad4db65SToby Isaac 4574165533cSJose E. Roman Input Parameters: 458fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 459fad4db65SToby Isaac . k - an integer in [0, n] 460fad4db65SToby Isaac - j - an index in [0, n choose k) 461fad4db65SToby Isaac 4624165533cSJose E. Roman Output Parameters: 463fad4db65SToby Isaac + perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set. 46428222859SToby Isaac - isOdd - if not NULL, return whether perm is an even or odd permutation. 465fad4db65SToby Isaac 466*87497f52SBarry Smith Note: 467*87497f52SBarry Smith Limited by arguments such that n choose k can be represented by `PetscInt` 468fad4db65SToby Isaac 469fad4db65SToby Isaac Level: beginner 470fad4db65SToby Isaac 471db781477SPatrick Sanan .seealso: `PetscDTEnumSubset()`, `PetscDTSubsetIndex()` 472fad4db65SToby Isaac M*/ 4739fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd) 4741a989b97SToby Isaac { 4755f80ce2aSJacob Faibussowitsch PetscInt i, l, m, Nk, odd = 0; 4765f80ce2aSJacob Faibussowitsch PetscInt *subcomp = perm+k; 4771a989b97SToby Isaac 47828222859SToby Isaac PetscFunctionBegin; 47928222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 4809566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 4811a989b97SToby Isaac for (i = 0, l = 0, m = 0; i < n && l < k; i++) { 4821a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 4831a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 4841a989b97SToby Isaac 4851a989b97SToby Isaac if (j < Nminuskminus) { 486fad4db65SToby Isaac perm[l++] = i; 4871a989b97SToby Isaac Nk = Nminuskminus; 4881a989b97SToby Isaac } else { 4891a989b97SToby Isaac subcomp[m++] = i; 4901a989b97SToby Isaac j -= Nminuskminus; 4911a989b97SToby Isaac odd ^= ((k - l) & 1); 4921a989b97SToby Isaac Nk = Nminusk; 4931a989b97SToby Isaac } 4941a989b97SToby Isaac } 4955f80ce2aSJacob Faibussowitsch for (; i < n; i++) subcomp[m++] = i; 4961a989b97SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 4971a989b97SToby Isaac PetscFunctionReturn(0); 4981a989b97SToby Isaac } 4991a989b97SToby Isaac 500ef0bb6c7SMatthew G. Knepley struct _p_PetscTabulation { 501a5b23f4aSJose E. Roman PetscInt K; /* Indicates a k-jet, namely tabulated derivatives up to order k */ 50219815104SMartin Diehl PetscInt Nr; /* The number of tabulation replicas (often 1) */ 503ef0bb6c7SMatthew G. Knepley PetscInt Np; /* The number of tabulation points in a replica */ 504ef0bb6c7SMatthew G. Knepley PetscInt Nb; /* The number of functions tabulated */ 505ef0bb6c7SMatthew G. Knepley PetscInt Nc; /* The number of function components */ 506ef0bb6c7SMatthew G. Knepley PetscInt cdim; /* The coordinate dimension */ 507ef0bb6c7SMatthew G. Knepley PetscReal **T; /* The tabulation T[K] of functions and their derivatives 508ef0bb6c7SMatthew G. Knepley T[0] = B[Nr*Np][Nb][Nc]: The basis function values at quadrature points 509ef0bb6c7SMatthew G. Knepley T[1] = D[Nr*Np][Nb][Nc][cdim]: The basis function derivatives at quadrature points 510ef0bb6c7SMatthew G. Knepley T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */ 511ef0bb6c7SMatthew G. Knepley }; 512ef0bb6c7SMatthew G. Knepley typedef struct _p_PetscTabulation *PetscTabulation; 513ef0bb6c7SMatthew G. Knepley 514d6685f55SMatthew G. Knepley typedef PetscErrorCode (*PetscProbFunc)(const PetscReal[], const PetscReal[], PetscReal[]); 515d6685f55SMatthew G. Knepley 516d6685f55SMatthew G. Knepley typedef enum {DTPROB_DENSITY_CONSTANT, DTPROB_DENSITY_GAUSSIAN, DTPROB_DENSITY_MAXWELL_BOLTZMANN, DTPROB_NUM_DENSITY} DTProbDensityType; 517d6685f55SMatthew G. Knepley PETSC_EXTERN const char * const DTProbDensityTypes[]; 518d6685f55SMatthew G. Knepley 519d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]); 520d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]); 521d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]); 522d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]); 523d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]); 524d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]); 525d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 526d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 527d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 528d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]); 529d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]); 530ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]); 531ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]); 532d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 533d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 534d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 535ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 536ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 537ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 538ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 539ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 540ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 541d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFunc *, PetscProbFunc *, PetscProbFunc *); 542d6685f55SMatthew G. Knepley 543d6685f55SMatthew G. Knepley #include <petscvec.h> 544d6685f55SMatthew G. Knepley 545d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFunc, PetscReal *); 546d6685f55SMatthew G. Knepley 54737045ce4SJed Brown #endif 548