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 92cd22861SMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCQUADRATURE_CLASSID; 102cd22861SMatthew G. Knepley 1121454ff5SMatthew G. Knepley /*S 1221454ff5SMatthew G. Knepley PetscQuadrature - Quadrature rule for integration. 1321454ff5SMatthew G. Knepley 14329bbf4eSMatthew G. Knepley Level: beginner 1521454ff5SMatthew G. Knepley 1621454ff5SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy() 1721454ff5SMatthew G. Knepley S*/ 1821454ff5SMatthew G. Knepley typedef struct _p_PetscQuadrature *PetscQuadrature; 1921454ff5SMatthew G. Knepley 208272889dSSatish Balay /*E 21916e780bShannah_mairs PetscGaussLobattoLegendreCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights 228272889dSSatish Balay 238272889dSSatish Balay Level: intermediate 248272889dSSatish Balay 25f2e8fe4dShannah_mairs $ PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA - compute the nodes via linear algebra 26d410ae54Shannah_mairs $ PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON - compute the nodes by solving a nonlinear equation with Newton's method 278272889dSSatish Balay 288272889dSSatish Balay E*/ 29f2e8fe4dShannah_mairs typedef enum {PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON} PetscGaussLobattoLegendreCreateType; 308272889dSSatish Balay 31d4afb720SToby Isaac /*E 32d4afb720SToby Isaac PetscDTNodeType - A description of strategies for generating nodes (both 33d4afb720SToby Isaac quadrature nodes and nodes for Lagrange polynomials) 34d4afb720SToby Isaac 35d4afb720SToby Isaac Level: intermediate 36d4afb720SToby Isaac 37d4afb720SToby Isaac $ PETSCDTNODES_DEFAULT - Nodes chosen by PETSc 38d4afb720SToby Isaac $ PETSCDTNODES_GAUSSJACOBI - Nodes at either Gauss-Jacobi or Gauss-Lobatto-Jacobi quadrature points 39d4afb720SToby Isaac $ PETSCDTNODES_EQUISPACED - Nodes equispaced either including the endpoints or excluding them 40d4afb720SToby Isaac $ PETSCDTNODES_TANHSINH - Nodes at Tanh-Sinh quadrature points 41d4afb720SToby Isaac 42d4afb720SToby Isaac Note: a PetscDTNodeType can be paired with a PetscBool to indicate whether 43d4afb720SToby Isaac the nodes include endpoints or not, and in the case of PETSCDT_GAUSSJACOBI 44d4afb720SToby Isaac with exponents for the weight function. 45d4afb720SToby Isaac 46d4afb720SToby Isaac E*/ 47d4afb720SToby Isaac typedef enum {PETSCDTNODES_DEFAULT=-1, PETSCDTNODES_GAUSSJACOBI, PETSCDTNODES_EQUISPACED, PETSCDTNODES_TANHSINH} PetscDTNodeType; 48d4afb720SToby Isaac 49d4afb720SToby Isaac PETSC_EXTERN const char *const PetscDTNodeTypes[]; 50d4afb720SToby Isaac 5121454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *); 52c9638911SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *); 53bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt*); 54bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt); 55a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt*); 56a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt); 574f9ab2b4SJed Brown PETSC_EXTERN PetscErrorCode PetscQuadratureEqual(PetscQuadrature, PetscQuadrature, PetscBool*); 58a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt*, PetscInt*, PetscInt*, const PetscReal *[], const PetscReal *[]); 59a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal [], const PetscReal []); 6021454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer); 6121454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *); 62a0845e3aSMatthew G. Knepley 632df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTensorQuadratureCreate(PetscQuadrature, PetscQuadrature, PetscQuadrature *); 6489710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *); 6589710940SMatthew G. Knepley 66907761f8SToby Isaac PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *); 67907761f8SToby Isaac 6837045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*); 69fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiNorm(PetscReal,PetscReal,PetscInt,PetscReal *); 7094e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt,PetscReal,PetscReal,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*); 71fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEvalJet(PetscReal,PetscReal,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscReal[]); 72fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPKDEvalJet(PetscInt,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscReal[]); 73d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedSize(PetscInt,PetscInt,PetscInt,PetscInt*); 74d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscInt,PetscReal[]); 7537045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt,PetscReal,PetscReal,PetscReal*,PetscReal*); 7694e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*,PetscReal*); 7794e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*,PetscReal*); 78916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt,PetscGaussLobattoLegendreCreateType,PetscReal*,PetscReal*); 79194825f6SJed Brown PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt,PetscInt,const PetscReal*,PetscInt,const PetscReal*,PetscReal*); 80a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*); 81e6a796c3SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*); 8237045ce4SJed Brown 83b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *); 84d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *); 85d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *); 86b3c0f97bSTom Klotz 87916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *); 88916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 89916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 90916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 91916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 92916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 93916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 94916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 95916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 96916e780bShannah_mairs 971a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 981a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 991a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 1001a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *); 1011a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *); 1021a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 1031a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *); 104dda711d0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]); 1051a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 1061a989b97SToby Isaac 107d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt,PetscInt,const PetscInt[],PetscInt*); 108d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt,PetscInt,PetscInt,PetscInt[]); 109fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGradedOrderToIndex(PetscInt,const PetscInt[],PetscInt*); 110fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToGradedOrder(PetscInt,PetscInt,PetscInt[]); 111d4afb720SToby Isaac 112fad4db65SToby Isaac #if defined(PETSC_USE_64BIT_INDICES) 113fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 20 114fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX 61 115fad4db65SToby Isaac #else 116fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 12 117fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX 29 118fad4db65SToby Isaac #endif 119fad4db65SToby Isaac 120fad4db65SToby Isaac /*MC 121fad4db65SToby Isaac PetscDTFactorial - Approximate n! as a real number 122fad4db65SToby Isaac 1234165533cSJose E. Roman Input Parameter: 124fad4db65SToby Isaac . n - a non-negative integer 125fad4db65SToby Isaac 1264165533cSJose E. Roman Output Parameter: 127fad4db65SToby Isaac . factorial - n! 128fad4db65SToby Isaac 129fad4db65SToby Isaac Level: beginner 130fad4db65SToby Isaac M*/ 1319fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial) 132fad4db65SToby Isaac { 133fad4db65SToby Isaac PetscReal f = 1.0; 134fad4db65SToby Isaac 135fad4db65SToby Isaac PetscFunctionBegin; 136e2ab39ccSLisandro Dalcin *factorial = -1.0; 1372c71b3e2SJacob Faibussowitsch PetscCheck(n >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %D", n); 1385f80ce2aSJacob Faibussowitsch for (PetscInt i = 1; i < n+1; ++i) f *= (PetscReal)i; 139fad4db65SToby Isaac *factorial = f; 140fad4db65SToby Isaac PetscFunctionReturn(0); 141fad4db65SToby Isaac } 142fad4db65SToby Isaac 143fad4db65SToby Isaac /*MC 144fad4db65SToby Isaac PetscDTFactorialInt - Compute n! as an integer 145fad4db65SToby Isaac 1464165533cSJose E. Roman Input Parameter: 147fad4db65SToby Isaac . n - a non-negative integer 148fad4db65SToby Isaac 1494165533cSJose E. Roman Output Parameter: 150fad4db65SToby Isaac . factorial - n! 151fad4db65SToby Isaac 152fad4db65SToby Isaac Level: beginner 153fad4db65SToby Isaac 154fad4db65SToby 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. 155fad4db65SToby Isaac M*/ 1569fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial) 157fad4db65SToby Isaac { 158fad4db65SToby Isaac PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600}; 159fad4db65SToby Isaac 16028222859SToby Isaac PetscFunctionBegin; 16128222859SToby Isaac *factorial = -1; 1622c71b3e2SJacob Faibussowitsch PetscCheck(n >= 0 && n <= PETSC_FACTORIAL_MAX,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of elements %D is not in supported range [0,%D]",n,PETSC_FACTORIAL_MAX); 163fad4db65SToby Isaac if (n <= 12) { 164fad4db65SToby Isaac *factorial = facLookup[n]; 165fad4db65SToby Isaac } else { 166fad4db65SToby Isaac PetscInt f = facLookup[12]; 167fad4db65SToby Isaac PetscInt i; 168fad4db65SToby Isaac 169fad4db65SToby Isaac for (i = 13; i < n+1; ++i) f *= i; 170fad4db65SToby Isaac *factorial = f; 171fad4db65SToby Isaac } 172fad4db65SToby Isaac PetscFunctionReturn(0); 173fad4db65SToby Isaac } 174fad4db65SToby Isaac 175fad4db65SToby Isaac /*MC 176fad4db65SToby Isaac PetscDTBinomial - Approximate the binomial coefficient "n choose k" 177fad4db65SToby Isaac 1784165533cSJose E. Roman Input Parameters: 179fad4db65SToby Isaac + n - a non-negative integer 180fad4db65SToby Isaac - k - an integer between 0 and n, inclusive 181fad4db65SToby Isaac 1824165533cSJose E. Roman Output Parameter: 183fad4db65SToby Isaac . binomial - approximation of the binomial coefficient n choose k 184fad4db65SToby Isaac 185fad4db65SToby Isaac Level: beginner 186fad4db65SToby Isaac M*/ 1879fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial) 1881a989b97SToby Isaac { 1891a989b97SToby Isaac PetscFunctionBeginHot; 190e2ab39ccSLisandro Dalcin *binomial = -1.0; 1912c71b3e2SJacob Faibussowitsch PetscCheck(n >= 0 && k >= 0 && k <= n,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial arguments (%D %D) must be non-negative, k <= n", n, k); 1921a989b97SToby Isaac if (n <= 3) { 1931a989b97SToby Isaac PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}}; 1941a989b97SToby Isaac 195e2ab39ccSLisandro Dalcin *binomial = (PetscReal)binomLookup[n][k]; 1961a989b97SToby Isaac } else { 197e2ab39ccSLisandro Dalcin PetscReal binom = 1.0; 1981a989b97SToby Isaac 1991a989b97SToby Isaac k = PetscMin(k, n - k); 2005f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1); 2011a989b97SToby Isaac *binomial = binom; 2021a989b97SToby Isaac } 2031a989b97SToby Isaac PetscFunctionReturn(0); 2041a989b97SToby Isaac } 2051a989b97SToby Isaac 206fad4db65SToby Isaac /*MC 207fad4db65SToby Isaac PetscDTBinomialInt - Compute the binomial coefficient "n choose k" 208fad4db65SToby Isaac 20997bb3fdcSJose E. Roman Input Parameters: 210fad4db65SToby Isaac + n - a non-negative integer 211fad4db65SToby Isaac - k - an integer between 0 and n, inclusive 212fad4db65SToby Isaac 21397bb3fdcSJose E. Roman Output Parameter: 214fad4db65SToby Isaac . binomial - the binomial coefficient n choose k 215fad4db65SToby Isaac 216fad4db65SToby Isaac Note: this is limited by integers that can be represented by PetscInt 217fad4db65SToby Isaac 218fad4db65SToby Isaac Level: beginner 219fad4db65SToby Isaac M*/ 2209fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial) 221fad4db65SToby Isaac { 22228222859SToby Isaac PetscInt bin; 22328222859SToby Isaac 22428222859SToby Isaac PetscFunctionBegin; 22528222859SToby Isaac *binomial = -1; 2262c71b3e2SJacob Faibussowitsch PetscCheck(n >= 0 && k >= 0 && k <= n,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial arguments (%D %D) must be non-negative, k <= n", n, k); 2272c71b3e2SJacob Faibussowitsch PetscCheck(n <= PETSC_BINOMIAL_MAX,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial elements %D is larger than max for PetscInt, %D", n, PETSC_BINOMIAL_MAX); 228fad4db65SToby Isaac if (n <= 3) { 229fad4db65SToby Isaac PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}}; 230fad4db65SToby Isaac 23128222859SToby Isaac bin = binomLookup[n][k]; 232fad4db65SToby Isaac } else { 233fad4db65SToby Isaac PetscInt binom = 1; 234fad4db65SToby Isaac 235fad4db65SToby Isaac k = PetscMin(k, n - k); 2365f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1); 23728222859SToby Isaac bin = binom; 238fad4db65SToby Isaac } 23928222859SToby Isaac *binomial = bin; 240fad4db65SToby Isaac PetscFunctionReturn(0); 241fad4db65SToby Isaac } 242fad4db65SToby Isaac 243fad4db65SToby Isaac /*MC 244fad4db65SToby Isaac PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps. 245fad4db65SToby Isaac 246fad4db65SToby Isaac A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation, 247fad4db65SToby Isaac by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in 24828222859SToby Isaac some position j >= i. This swap is encoded as the difference (j - i). The difference d_i at step i is less than 24928222859SToby Isaac (n - i). This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number 250fad4db65SToby Isaac (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}. 251fad4db65SToby Isaac 2524165533cSJose E. Roman Input Parameters: 253fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 2548cd1e013SToby Isaac - k - an integer in [0, n!) 255fad4db65SToby Isaac 2564165533cSJose E. Roman Output Parameters: 257fad4db65SToby Isaac + perm - the permuted list of the integers [0, ..., n-1] 2588cd1e013SToby Isaac - isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps. 259fad4db65SToby Isaac 260fad4db65SToby 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. 261fad4db65SToby Isaac 262fad4db65SToby Isaac Level: beginner 263fad4db65SToby Isaac M*/ 2649fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd) 2651a989b97SToby Isaac { 2661a989b97SToby Isaac PetscInt odd = 0; 2671a989b97SToby Isaac PetscInt i; 268fad4db65SToby Isaac PetscInt work[PETSC_FACTORIAL_MAX]; 269fad4db65SToby Isaac PetscInt *w; 2701a989b97SToby Isaac 27128222859SToby Isaac PetscFunctionBegin; 27228222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 2732c71b3e2SJacob Faibussowitsch PetscCheck(n >= 0 && n <= PETSC_FACTORIAL_MAX,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of elements %D is not in supported range [0,%D]",n,PETSC_FACTORIAL_MAX); 274fad4db65SToby Isaac w = &work[n - 2]; 2751a989b97SToby Isaac for (i = 2; i <= n; i++) { 2761a989b97SToby Isaac *(w--) = k % i; 2771a989b97SToby Isaac k /= i; 2781a989b97SToby Isaac } 2791a989b97SToby Isaac for (i = 0; i < n; i++) perm[i] = i; 2801a989b97SToby Isaac for (i = 0; i < n - 1; i++) { 2811a989b97SToby Isaac PetscInt s = work[i]; 2821a989b97SToby Isaac PetscInt swap = perm[i]; 2831a989b97SToby Isaac 2841a989b97SToby Isaac perm[i] = perm[i + s]; 2851a989b97SToby Isaac perm[i + s] = swap; 2861a989b97SToby Isaac odd ^= (!!s); 2871a989b97SToby Isaac } 2881a989b97SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 2891a989b97SToby Isaac PetscFunctionReturn(0); 2901a989b97SToby Isaac } 2911a989b97SToby Isaac 292fad4db65SToby Isaac /*MC 2938cd1e013SToby Isaac PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!). This inverts PetscDTEnumPerm. 2948cd1e013SToby Isaac 2954165533cSJose E. Roman Input Parameters: 2968cd1e013SToby Isaac + n - a non-negative integer (see note about limits below) 2978cd1e013SToby Isaac - perm - the permuted list of the integers [0, ..., n-1] 2988cd1e013SToby Isaac 2994165533cSJose E. Roman Output Parameters: 3008cd1e013SToby Isaac + k - an integer in [0, n!) 301f0fc11ceSJed Brown - isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps. 3028cd1e013SToby Isaac 3038cd1e013SToby 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. 3048cd1e013SToby Isaac 3058cd1e013SToby Isaac Level: beginner 3068cd1e013SToby Isaac M*/ 3079fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd) 3088cd1e013SToby Isaac { 3098cd1e013SToby Isaac PetscInt odd = 0; 3108cd1e013SToby Isaac PetscInt i, idx; 3118cd1e013SToby Isaac PetscInt work[PETSC_FACTORIAL_MAX]; 3128cd1e013SToby Isaac PetscInt iwork[PETSC_FACTORIAL_MAX]; 3138cd1e013SToby Isaac 3148cd1e013SToby Isaac PetscFunctionBeginHot; 31528222859SToby Isaac *k = -1; 31628222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 3172c71b3e2SJacob Faibussowitsch PetscCheck(n >= 0 && n <= PETSC_FACTORIAL_MAX,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of elements %D is not in supported range [0,%D]",n,PETSC_FACTORIAL_MAX); 3188cd1e013SToby Isaac for (i = 0; i < n; i++) work[i] = i; /* partial permutation */ 3198cd1e013SToby Isaac for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */ 3208cd1e013SToby Isaac for (idx = 0, i = 0; i < n - 1; i++) { 3218cd1e013SToby Isaac PetscInt j = perm[i]; 3228cd1e013SToby Isaac PetscInt icur = work[i]; 3238cd1e013SToby Isaac PetscInt jloc = iwork[j]; 3248cd1e013SToby Isaac PetscInt diff = jloc - i; 3258cd1e013SToby Isaac 3268cd1e013SToby Isaac idx = idx * (n - i) + diff; 3278cd1e013SToby Isaac /* swap (i, jloc) */ 3288cd1e013SToby Isaac work[i] = j; 3298cd1e013SToby Isaac work[jloc] = icur; 3308cd1e013SToby Isaac iwork[j] = i; 3318cd1e013SToby Isaac iwork[icur] = jloc; 3328cd1e013SToby Isaac odd ^= (!!diff); 3338cd1e013SToby Isaac } 3348cd1e013SToby Isaac *k = idx; 3358cd1e013SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 3368cd1e013SToby Isaac PetscFunctionReturn(0); 3378cd1e013SToby Isaac } 3388cd1e013SToby Isaac 3398cd1e013SToby Isaac /*MC 340fad4db65SToby Isaac PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k). 341fad4db65SToby Isaac The encoding is in lexicographic order. 342fad4db65SToby Isaac 3434165533cSJose E. Roman Input Parameters: 344fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 345fad4db65SToby Isaac . k - an integer in [0, n] 346fad4db65SToby Isaac - j - an index in [0, n choose k) 347fad4db65SToby Isaac 3484165533cSJose E. Roman Output Parameter: 349fad4db65SToby Isaac . subset - the jth subset of size k of the integers [0, ..., n - 1] 350fad4db65SToby Isaac 351fad4db65SToby Isaac Note: this is limited by arguments such that n choose k can be represented by PetscInt 352fad4db65SToby Isaac 353fad4db65SToby Isaac Level: beginner 354fad4db65SToby Isaac 355fad4db65SToby Isaac .seealso: PetscDTSubsetIndex() 356fad4db65SToby Isaac M*/ 3579fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset) 3581a989b97SToby Isaac { 3595f80ce2aSJacob Faibussowitsch PetscInt Nk; 3601a989b97SToby Isaac 3611a989b97SToby Isaac PetscFunctionBeginHot; 3629566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 3635f80ce2aSJacob Faibussowitsch for (PetscInt i = 0, l = 0; i < n && l < k; i++) { 3641a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 3651a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 3661a989b97SToby Isaac 3671a989b97SToby Isaac if (j < Nminuskminus) { 3681a989b97SToby Isaac subset[l++] = i; 3691a989b97SToby Isaac Nk = Nminuskminus; 3701a989b97SToby Isaac } else { 3711a989b97SToby Isaac j -= Nminuskminus; 3721a989b97SToby Isaac Nk = Nminusk; 3731a989b97SToby Isaac } 3741a989b97SToby Isaac } 3751a989b97SToby Isaac PetscFunctionReturn(0); 3761a989b97SToby Isaac } 3771a989b97SToby Isaac 378fad4db65SToby Isaac /*MC 379fad4db65SToby Isaac 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. This is the inverse of PetscDTEnumSubset. 380fad4db65SToby Isaac 3814165533cSJose E. Roman Input Parameters: 382fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 383fad4db65SToby Isaac . k - an integer in [0, n] 384fad4db65SToby Isaac - subset - an ordered subset of the integers [0, ..., n - 1] 385fad4db65SToby Isaac 3864165533cSJose E. Roman Output Parameter: 387fad4db65SToby Isaac . index - the rank of the subset in lexicographic order 388fad4db65SToby Isaac 389fad4db65SToby Isaac Note: this is limited by arguments such that n choose k can be represented by PetscInt 390fad4db65SToby Isaac 391fad4db65SToby Isaac Level: beginner 392fad4db65SToby Isaac 393fad4db65SToby Isaac .seealso: PetscDTEnumSubset() 394fad4db65SToby Isaac M*/ 3959fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index) 3961a989b97SToby Isaac { 3975f80ce2aSJacob Faibussowitsch PetscInt j = 0, Nk; 3981a989b97SToby Isaac 39928222859SToby Isaac PetscFunctionBegin; 40028222859SToby Isaac *index = -1; 4019566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 4025f80ce2aSJacob Faibussowitsch for (PetscInt i = 0, l = 0; i < n && l < k; i++) { 4031a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 4041a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 4051a989b97SToby Isaac 4061a989b97SToby Isaac if (subset[l] == i) { 4071a989b97SToby Isaac l++; 4081a989b97SToby Isaac Nk = Nminuskminus; 4091a989b97SToby Isaac } else { 4101a989b97SToby Isaac j += Nminuskminus; 4111a989b97SToby Isaac Nk = Nminusk; 4121a989b97SToby Isaac } 4131a989b97SToby Isaac } 4141a989b97SToby Isaac *index = j; 4151a989b97SToby Isaac PetscFunctionReturn(0); 4161a989b97SToby Isaac } 4171a989b97SToby Isaac 418fad4db65SToby Isaac /*MC 41928222859SToby 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. 420fad4db65SToby Isaac 4214165533cSJose E. Roman Input Parameters: 422fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 423fad4db65SToby Isaac . k - an integer in [0, n] 424fad4db65SToby Isaac - j - an index in [0, n choose k) 425fad4db65SToby Isaac 4264165533cSJose E. Roman Output Parameters: 427fad4db65SToby Isaac + perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set. 42828222859SToby Isaac - isOdd - if not NULL, return whether perm is an even or odd permutation. 429fad4db65SToby Isaac 430fad4db65SToby Isaac Note: this is limited by arguments such that n choose k can be represented by PetscInt 431fad4db65SToby Isaac 432fad4db65SToby Isaac Level: beginner 433fad4db65SToby Isaac 434fad4db65SToby Isaac .seealso: PetscDTEnumSubset(), PetscDTSubsetIndex() 435fad4db65SToby Isaac M*/ 4369fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd) 4371a989b97SToby Isaac { 4385f80ce2aSJacob Faibussowitsch PetscInt i, l, m, Nk, odd = 0; 4395f80ce2aSJacob Faibussowitsch PetscInt *subcomp = perm+k; 4401a989b97SToby Isaac 44128222859SToby Isaac PetscFunctionBegin; 44228222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 4439566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n, k, &Nk)); 4441a989b97SToby Isaac for (i = 0, l = 0, m = 0; i < n && l < k; i++) { 4451a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 4461a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 4471a989b97SToby Isaac 4481a989b97SToby Isaac if (j < Nminuskminus) { 449fad4db65SToby Isaac perm[l++] = i; 4501a989b97SToby Isaac Nk = Nminuskminus; 4511a989b97SToby Isaac } else { 4521a989b97SToby Isaac subcomp[m++] = i; 4531a989b97SToby Isaac j -= Nminuskminus; 4541a989b97SToby Isaac odd ^= ((k - l) & 1); 4551a989b97SToby Isaac Nk = Nminusk; 4561a989b97SToby Isaac } 4571a989b97SToby Isaac } 4585f80ce2aSJacob Faibussowitsch for (; i < n; i++) subcomp[m++] = i; 4591a989b97SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 4601a989b97SToby Isaac PetscFunctionReturn(0); 4611a989b97SToby Isaac } 4621a989b97SToby Isaac 463ef0bb6c7SMatthew G. Knepley struct _p_PetscTabulation { 464a5b23f4aSJose E. Roman PetscInt K; /* Indicates a k-jet, namely tabulated derivatives up to order k */ 46519815104SMartin Diehl PetscInt Nr; /* The number of tabulation replicas (often 1) */ 466ef0bb6c7SMatthew G. Knepley PetscInt Np; /* The number of tabulation points in a replica */ 467ef0bb6c7SMatthew G. Knepley PetscInt Nb; /* The number of functions tabulated */ 468ef0bb6c7SMatthew G. Knepley PetscInt Nc; /* The number of function components */ 469ef0bb6c7SMatthew G. Knepley PetscInt cdim; /* The coordinate dimension */ 470ef0bb6c7SMatthew G. Knepley PetscReal **T; /* The tabulation T[K] of functions and their derivatives 471ef0bb6c7SMatthew G. Knepley T[0] = B[Nr*Np][Nb][Nc]: The basis function values at quadrature points 472ef0bb6c7SMatthew G. Knepley T[1] = D[Nr*Np][Nb][Nc][cdim]: The basis function derivatives at quadrature points 473ef0bb6c7SMatthew G. Knepley T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */ 474ef0bb6c7SMatthew G. Knepley }; 475ef0bb6c7SMatthew G. Knepley typedef struct _p_PetscTabulation *PetscTabulation; 476ef0bb6c7SMatthew G. Knepley 477d6685f55SMatthew G. Knepley typedef PetscErrorCode (*PetscProbFunc)(const PetscReal[], const PetscReal[], PetscReal[]); 478d6685f55SMatthew G. Knepley 479d6685f55SMatthew G. Knepley typedef enum {DTPROB_DENSITY_CONSTANT, DTPROB_DENSITY_GAUSSIAN, DTPROB_DENSITY_MAXWELL_BOLTZMANN, DTPROB_NUM_DENSITY} DTProbDensityType; 480d6685f55SMatthew G. Knepley PETSC_EXTERN const char * const DTProbDensityTypes[]; 481d6685f55SMatthew G. Knepley 482d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]); 483d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]); 484d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]); 485d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]); 486d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]); 487d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]); 488d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 489d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 490d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 491d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]); 492d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]); 493*ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]); 494*ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]); 495d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 496d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 497d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 498*ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 499*ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 500*ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant2D(const PetscReal[], const PetscReal[], PetscReal[]); 501*ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 502*ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 503*ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant3D(const PetscReal[], const PetscReal[], PetscReal[]); 504d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFunc *, PetscProbFunc *, PetscProbFunc *); 505d6685f55SMatthew G. Knepley 506d6685f55SMatthew G. Knepley #include <petscvec.h> 507d6685f55SMatthew G. Knepley 508d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFunc, PetscReal *); 509d6685f55SMatthew G. Knepley 51037045ce4SJed Brown #endif 511