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); 57*4f9ab2b4SJed 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 PetscInt i; 135fad4db65SToby Isaac 136fad4db65SToby Isaac PetscFunctionBegin; 137e2ab39ccSLisandro Dalcin *factorial = -1.0; 1382c71b3e2SJacob Faibussowitsch PetscCheck(n >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %D", n); 139e2ab39ccSLisandro Dalcin for (i = 1; i < n+1; ++i) f *= (PetscReal)i; 140fad4db65SToby Isaac *factorial = f; 141fad4db65SToby Isaac PetscFunctionReturn(0); 142fad4db65SToby Isaac } 143fad4db65SToby Isaac 144fad4db65SToby Isaac /*MC 145fad4db65SToby Isaac PetscDTFactorialInt - Compute n! as an integer 146fad4db65SToby Isaac 1474165533cSJose E. Roman Input Parameter: 148fad4db65SToby Isaac . n - a non-negative integer 149fad4db65SToby Isaac 1504165533cSJose E. Roman Output Parameter: 151fad4db65SToby Isaac . factorial - n! 152fad4db65SToby Isaac 153fad4db65SToby Isaac Level: beginner 154fad4db65SToby Isaac 155fad4db65SToby 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. 156fad4db65SToby Isaac M*/ 1579fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial) 158fad4db65SToby Isaac { 159fad4db65SToby Isaac PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600}; 160fad4db65SToby Isaac 16128222859SToby Isaac PetscFunctionBegin; 16228222859SToby Isaac *factorial = -1; 1632c71b3e2SJacob 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); 164fad4db65SToby Isaac if (n <= 12) { 165fad4db65SToby Isaac *factorial = facLookup[n]; 166fad4db65SToby Isaac } else { 167fad4db65SToby Isaac PetscInt f = facLookup[12]; 168fad4db65SToby Isaac PetscInt i; 169fad4db65SToby Isaac 170fad4db65SToby Isaac for (i = 13; i < n+1; ++i) f *= i; 171fad4db65SToby Isaac *factorial = f; 172fad4db65SToby Isaac } 173fad4db65SToby Isaac PetscFunctionReturn(0); 174fad4db65SToby Isaac } 175fad4db65SToby Isaac 176fad4db65SToby Isaac /*MC 177fad4db65SToby Isaac PetscDTBinomial - Approximate the binomial coefficient "n choose k" 178fad4db65SToby Isaac 1794165533cSJose E. Roman Input Parameters: 180fad4db65SToby Isaac + n - a non-negative integer 181fad4db65SToby Isaac - k - an integer between 0 and n, inclusive 182fad4db65SToby Isaac 1834165533cSJose E. Roman Output Parameter: 184fad4db65SToby Isaac . binomial - approximation of the binomial coefficient n choose k 185fad4db65SToby Isaac 186fad4db65SToby Isaac Level: beginner 187fad4db65SToby Isaac M*/ 1889fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial) 1891a989b97SToby Isaac { 1901a989b97SToby Isaac PetscFunctionBeginHot; 191e2ab39ccSLisandro Dalcin *binomial = -1.0; 1922c71b3e2SJacob 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); 1931a989b97SToby Isaac if (n <= 3) { 1941a989b97SToby Isaac PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}}; 1951a989b97SToby Isaac 196e2ab39ccSLisandro Dalcin *binomial = (PetscReal)binomLookup[n][k]; 1971a989b97SToby Isaac } else { 198e2ab39ccSLisandro Dalcin PetscReal binom = 1.0; 1991a989b97SToby Isaac PetscInt i; 2001a989b97SToby Isaac 2011a989b97SToby Isaac k = PetscMin(k, n - k); 202e2ab39ccSLisandro Dalcin for (i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1); 2031a989b97SToby Isaac *binomial = binom; 2041a989b97SToby Isaac } 2051a989b97SToby Isaac PetscFunctionReturn(0); 2061a989b97SToby Isaac } 2071a989b97SToby Isaac 208fad4db65SToby Isaac /*MC 209fad4db65SToby Isaac PetscDTBinomialInt - Compute the binomial coefficient "n choose k" 210fad4db65SToby Isaac 21197bb3fdcSJose E. Roman Input Parameters: 212fad4db65SToby Isaac + n - a non-negative integer 213fad4db65SToby Isaac - k - an integer between 0 and n, inclusive 214fad4db65SToby Isaac 21597bb3fdcSJose E. Roman Output Parameter: 216fad4db65SToby Isaac . binomial - the binomial coefficient n choose k 217fad4db65SToby Isaac 218fad4db65SToby Isaac Note: this is limited by integers that can be represented by PetscInt 219fad4db65SToby Isaac 220fad4db65SToby Isaac Level: beginner 221fad4db65SToby Isaac M*/ 2229fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial) 223fad4db65SToby Isaac { 22428222859SToby Isaac PetscInt bin; 22528222859SToby Isaac 22628222859SToby Isaac PetscFunctionBegin; 22728222859SToby Isaac *binomial = -1; 2282c71b3e2SJacob 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); 2292c71b3e2SJacob 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); 230fad4db65SToby Isaac if (n <= 3) { 231fad4db65SToby Isaac PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}}; 232fad4db65SToby Isaac 23328222859SToby Isaac bin = binomLookup[n][k]; 234fad4db65SToby Isaac } else { 235fad4db65SToby Isaac PetscInt binom = 1; 236fad4db65SToby Isaac PetscInt i; 237fad4db65SToby Isaac 238fad4db65SToby Isaac k = PetscMin(k, n - k); 239fad4db65SToby Isaac for (i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1); 24028222859SToby Isaac bin = binom; 241fad4db65SToby Isaac } 24228222859SToby Isaac *binomial = bin; 243fad4db65SToby Isaac PetscFunctionReturn(0); 244fad4db65SToby Isaac } 245fad4db65SToby Isaac 246fad4db65SToby Isaac /*MC 247fad4db65SToby Isaac PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps. 248fad4db65SToby Isaac 249fad4db65SToby Isaac A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation, 250fad4db65SToby Isaac by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in 25128222859SToby Isaac some position j >= i. This swap is encoded as the difference (j - i). The difference d_i at step i is less than 25228222859SToby Isaac (n - i). This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number 253fad4db65SToby Isaac (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}. 254fad4db65SToby Isaac 2554165533cSJose E. Roman Input Parameters: 256fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 2578cd1e013SToby Isaac - k - an integer in [0, n!) 258fad4db65SToby Isaac 2594165533cSJose E. Roman Output Parameters: 260fad4db65SToby Isaac + perm - the permuted list of the integers [0, ..., n-1] 2618cd1e013SToby Isaac - isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps. 262fad4db65SToby Isaac 263fad4db65SToby 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. 264fad4db65SToby Isaac 265fad4db65SToby Isaac Level: beginner 266fad4db65SToby Isaac M*/ 2679fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd) 2681a989b97SToby Isaac { 2691a989b97SToby Isaac PetscInt odd = 0; 2701a989b97SToby Isaac PetscInt i; 271fad4db65SToby Isaac PetscInt work[PETSC_FACTORIAL_MAX]; 272fad4db65SToby Isaac PetscInt *w; 2731a989b97SToby Isaac 27428222859SToby Isaac PetscFunctionBegin; 27528222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 2762c71b3e2SJacob 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); 277fad4db65SToby Isaac w = &work[n - 2]; 2781a989b97SToby Isaac for (i = 2; i <= n; i++) { 2791a989b97SToby Isaac *(w--) = k % i; 2801a989b97SToby Isaac k /= i; 2811a989b97SToby Isaac } 2821a989b97SToby Isaac for (i = 0; i < n; i++) perm[i] = i; 2831a989b97SToby Isaac for (i = 0; i < n - 1; i++) { 2841a989b97SToby Isaac PetscInt s = work[i]; 2851a989b97SToby Isaac PetscInt swap = perm[i]; 2861a989b97SToby Isaac 2871a989b97SToby Isaac perm[i] = perm[i + s]; 2881a989b97SToby Isaac perm[i + s] = swap; 2891a989b97SToby Isaac odd ^= (!!s); 2901a989b97SToby Isaac } 2911a989b97SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 2921a989b97SToby Isaac PetscFunctionReturn(0); 2931a989b97SToby Isaac } 2941a989b97SToby Isaac 295fad4db65SToby Isaac /*MC 2968cd1e013SToby Isaac PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!). This inverts PetscDTEnumPerm. 2978cd1e013SToby Isaac 2984165533cSJose E. Roman Input Parameters: 2998cd1e013SToby Isaac + n - a non-negative integer (see note about limits below) 3008cd1e013SToby Isaac - perm - the permuted list of the integers [0, ..., n-1] 3018cd1e013SToby Isaac 3024165533cSJose E. Roman Output Parameters: 3038cd1e013SToby Isaac + k - an integer in [0, n!) 304f0fc11ceSJed Brown - isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps. 3058cd1e013SToby Isaac 3068cd1e013SToby 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. 3078cd1e013SToby Isaac 3088cd1e013SToby Isaac Level: beginner 3098cd1e013SToby Isaac M*/ 3109fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd) 3118cd1e013SToby Isaac { 3128cd1e013SToby Isaac PetscInt odd = 0; 3138cd1e013SToby Isaac PetscInt i, idx; 3148cd1e013SToby Isaac PetscInt work[PETSC_FACTORIAL_MAX]; 3158cd1e013SToby Isaac PetscInt iwork[PETSC_FACTORIAL_MAX]; 3168cd1e013SToby Isaac 3178cd1e013SToby Isaac PetscFunctionBeginHot; 31828222859SToby Isaac *k = -1; 31928222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 3202c71b3e2SJacob 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); 3218cd1e013SToby Isaac for (i = 0; i < n; i++) work[i] = i; /* partial permutation */ 3228cd1e013SToby Isaac for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */ 3238cd1e013SToby Isaac for (idx = 0, i = 0; i < n - 1; i++) { 3248cd1e013SToby Isaac PetscInt j = perm[i]; 3258cd1e013SToby Isaac PetscInt icur = work[i]; 3268cd1e013SToby Isaac PetscInt jloc = iwork[j]; 3278cd1e013SToby Isaac PetscInt diff = jloc - i; 3288cd1e013SToby Isaac 3298cd1e013SToby Isaac idx = idx * (n - i) + diff; 3308cd1e013SToby Isaac /* swap (i, jloc) */ 3318cd1e013SToby Isaac work[i] = j; 3328cd1e013SToby Isaac work[jloc] = icur; 3338cd1e013SToby Isaac iwork[j] = i; 3348cd1e013SToby Isaac iwork[icur] = jloc; 3358cd1e013SToby Isaac odd ^= (!!diff); 3368cd1e013SToby Isaac } 3378cd1e013SToby Isaac *k = idx; 3388cd1e013SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 3398cd1e013SToby Isaac PetscFunctionReturn(0); 3408cd1e013SToby Isaac } 3418cd1e013SToby Isaac 3428cd1e013SToby Isaac /*MC 343fad4db65SToby Isaac PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k). 344fad4db65SToby Isaac The encoding is in lexicographic order. 345fad4db65SToby Isaac 3464165533cSJose E. Roman Input Parameters: 347fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 348fad4db65SToby Isaac . k - an integer in [0, n] 349fad4db65SToby Isaac - j - an index in [0, n choose k) 350fad4db65SToby Isaac 3514165533cSJose E. Roman Output Parameter: 352fad4db65SToby Isaac . subset - the jth subset of size k of the integers [0, ..., n - 1] 353fad4db65SToby Isaac 354fad4db65SToby Isaac Note: this is limited by arguments such that n choose k can be represented by PetscInt 355fad4db65SToby Isaac 356fad4db65SToby Isaac Level: beginner 357fad4db65SToby Isaac 358fad4db65SToby Isaac .seealso: PetscDTSubsetIndex() 359fad4db65SToby Isaac M*/ 3609fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset) 3611a989b97SToby Isaac { 3621a989b97SToby Isaac PetscInt Nk, i, l; 3631a989b97SToby Isaac PetscErrorCode ierr; 3641a989b97SToby Isaac 3651a989b97SToby Isaac PetscFunctionBeginHot; 366fad4db65SToby Isaac ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr); 3671a989b97SToby Isaac for (i = 0, l = 0; i < n && l < k; i++) { 3681a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 3691a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 3701a989b97SToby Isaac 3711a989b97SToby Isaac if (j < Nminuskminus) { 3721a989b97SToby Isaac subset[l++] = i; 3731a989b97SToby Isaac Nk = Nminuskminus; 3741a989b97SToby Isaac } else { 3751a989b97SToby Isaac j -= Nminuskminus; 3761a989b97SToby Isaac Nk = Nminusk; 3771a989b97SToby Isaac } 3781a989b97SToby Isaac } 3791a989b97SToby Isaac PetscFunctionReturn(0); 3801a989b97SToby Isaac } 3811a989b97SToby Isaac 382fad4db65SToby Isaac /*MC 383fad4db65SToby 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. 384fad4db65SToby Isaac 3854165533cSJose E. Roman Input Parameters: 386fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 387fad4db65SToby Isaac . k - an integer in [0, n] 388fad4db65SToby Isaac - subset - an ordered subset of the integers [0, ..., n - 1] 389fad4db65SToby Isaac 3904165533cSJose E. Roman Output Parameter: 391fad4db65SToby Isaac . index - the rank of the subset in lexicographic order 392fad4db65SToby Isaac 393fad4db65SToby Isaac Note: this is limited by arguments such that n choose k can be represented by PetscInt 394fad4db65SToby Isaac 395fad4db65SToby Isaac Level: beginner 396fad4db65SToby Isaac 397fad4db65SToby Isaac .seealso: PetscDTEnumSubset() 398fad4db65SToby Isaac M*/ 3999fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index) 4001a989b97SToby Isaac { 4011a989b97SToby Isaac PetscInt i, j = 0, l, Nk; 4021a989b97SToby Isaac PetscErrorCode ierr; 4031a989b97SToby Isaac 40428222859SToby Isaac PetscFunctionBegin; 40528222859SToby Isaac *index = -1; 406fad4db65SToby Isaac ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr); 4071a989b97SToby Isaac for (i = 0, l = 0; i < n && l < k; i++) { 4081a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 4091a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 4101a989b97SToby Isaac 4111a989b97SToby Isaac if (subset[l] == i) { 4121a989b97SToby Isaac l++; 4131a989b97SToby Isaac Nk = Nminuskminus; 4141a989b97SToby Isaac } else { 4151a989b97SToby Isaac j += Nminuskminus; 4161a989b97SToby Isaac Nk = Nminusk; 4171a989b97SToby Isaac } 4181a989b97SToby Isaac } 4191a989b97SToby Isaac *index = j; 4201a989b97SToby Isaac PetscFunctionReturn(0); 4211a989b97SToby Isaac } 4221a989b97SToby Isaac 423fad4db65SToby Isaac /*MC 42428222859SToby 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. 425fad4db65SToby Isaac 4264165533cSJose E. Roman Input Parameters: 427fad4db65SToby Isaac + n - a non-negative integer (see note about limits below) 428fad4db65SToby Isaac . k - an integer in [0, n] 429fad4db65SToby Isaac - j - an index in [0, n choose k) 430fad4db65SToby Isaac 4314165533cSJose E. Roman Output Parameters: 432fad4db65SToby Isaac + perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set. 43328222859SToby Isaac - isOdd - if not NULL, return whether perm is an even or odd permutation. 434fad4db65SToby Isaac 435fad4db65SToby Isaac Note: this is limited by arguments such that n choose k can be represented by PetscInt 436fad4db65SToby Isaac 437fad4db65SToby Isaac Level: beginner 438fad4db65SToby Isaac 439fad4db65SToby Isaac .seealso: PetscDTEnumSubset(), PetscDTSubsetIndex() 440fad4db65SToby Isaac M*/ 4419fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd) 4421a989b97SToby Isaac { 4431a989b97SToby Isaac PetscInt i, l, m, *subcomp, Nk; 4441a989b97SToby Isaac PetscInt odd; 4451a989b97SToby Isaac PetscErrorCode ierr; 4461a989b97SToby Isaac 44728222859SToby Isaac PetscFunctionBegin; 44828222859SToby Isaac if (isOdd) *isOdd = PETSC_FALSE; 449fad4db65SToby Isaac ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr); 4501a989b97SToby Isaac odd = 0; 451fad4db65SToby Isaac subcomp = &perm[k]; 4521a989b97SToby Isaac for (i = 0, l = 0, m = 0; i < n && l < k; i++) { 4531a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 4541a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 4551a989b97SToby Isaac 4561a989b97SToby Isaac if (j < Nminuskminus) { 457fad4db65SToby Isaac perm[l++] = i; 4581a989b97SToby Isaac Nk = Nminuskminus; 4591a989b97SToby Isaac } else { 4601a989b97SToby Isaac subcomp[m++] = i; 4611a989b97SToby Isaac j -= Nminuskminus; 4621a989b97SToby Isaac odd ^= ((k - l) & 1); 4631a989b97SToby Isaac Nk = Nminusk; 4641a989b97SToby Isaac } 4651a989b97SToby Isaac } 4661a989b97SToby Isaac for (; i < n; i++) { 4671a989b97SToby Isaac subcomp[m++] = i; 4681a989b97SToby Isaac } 4691a989b97SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 4701a989b97SToby Isaac PetscFunctionReturn(0); 4711a989b97SToby Isaac } 4721a989b97SToby Isaac 473ef0bb6c7SMatthew G. Knepley struct _p_PetscTabulation { 474a5b23f4aSJose E. Roman PetscInt K; /* Indicates a k-jet, namely tabulated derivatives up to order k */ 47519815104SMartin Diehl PetscInt Nr; /* The number of tabulation replicas (often 1) */ 476ef0bb6c7SMatthew G. Knepley PetscInt Np; /* The number of tabulation points in a replica */ 477ef0bb6c7SMatthew G. Knepley PetscInt Nb; /* The number of functions tabulated */ 478ef0bb6c7SMatthew G. Knepley PetscInt Nc; /* The number of function components */ 479ef0bb6c7SMatthew G. Knepley PetscInt cdim; /* The coordinate dimension */ 480ef0bb6c7SMatthew G. Knepley PetscReal **T; /* The tabulation T[K] of functions and their derivatives 481ef0bb6c7SMatthew G. Knepley T[0] = B[Nr*Np][Nb][Nc]: The basis function values at quadrature points 482ef0bb6c7SMatthew G. Knepley T[1] = D[Nr*Np][Nb][Nc][cdim]: The basis function derivatives at quadrature points 483ef0bb6c7SMatthew G. Knepley T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */ 484ef0bb6c7SMatthew G. Knepley }; 485ef0bb6c7SMatthew G. Knepley typedef struct _p_PetscTabulation *PetscTabulation; 486ef0bb6c7SMatthew G. Knepley 487d6685f55SMatthew G. Knepley typedef PetscErrorCode (*PetscProbFunc)(const PetscReal[], const PetscReal[], PetscReal[]); 488d6685f55SMatthew G. Knepley 489d6685f55SMatthew G. Knepley typedef enum {DTPROB_DENSITY_CONSTANT, DTPROB_DENSITY_GAUSSIAN, DTPROB_DENSITY_MAXWELL_BOLTZMANN, DTPROB_NUM_DENSITY} DTProbDensityType; 490d6685f55SMatthew G. Knepley PETSC_EXTERN const char * const DTProbDensityTypes[]; 491d6685f55SMatthew G. Knepley 492d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]); 493d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]); 494d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]); 495d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]); 496d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]); 497d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]); 498d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 499d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 500d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]); 501d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]); 502d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]); 503d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 504d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 505d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]); 506d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFunc *, PetscProbFunc *, PetscProbFunc *); 507d6685f55SMatthew G. Knepley 508d6685f55SMatthew G. Knepley #include <petscvec.h> 509d6685f55SMatthew G. Knepley 510d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFunc, PetscReal *); 511d6685f55SMatthew G. Knepley 51237045ce4SJed Brown #endif 513