xref: /petsc/include/petscdt.h (revision ac09b9214d23ea9ad238aa607de9fa447fd4e91b)
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 
9*ac09b921SBarry Smith /* SUBMANSEC = DT */
10*ac09b921SBarry 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 
27f2e8fe4dShannah_mairs $  PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA - compute the nodes via linear algebra
28d410ae54Shannah_mairs $  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 
39d4afb720SToby Isaac $  PETSCDTNODES_DEFAULT - Nodes chosen by PETSc
40d4afb720SToby Isaac $  PETSCDTNODES_GAUSSJACOBI - Nodes at either Gauss-Jacobi or Gauss-Lobatto-Jacobi quadrature points
41d4afb720SToby Isaac $  PETSCDTNODES_EQUISPACED - Nodes equispaced either including the endpoints or excluding them
42d4afb720SToby Isaac $  PETSCDTNODES_TANHSINH - Nodes at Tanh-Sinh quadrature points
43d4afb720SToby Isaac 
44d4afb720SToby Isaac   Note: a PetscDTNodeType can be paired with a PetscBool to indicate whether
45d4afb720SToby Isaac   the nodes include endpoints or not, and in the case of PETSCDT_GAUSSJACOBI
46d4afb720SToby Isaac   with exponents for the weight function.
47d4afb720SToby Isaac 
48d4afb720SToby Isaac E*/
49d4afb720SToby Isaac typedef enum {PETSCDTNODES_DEFAULT=-1, PETSCDTNODES_GAUSSJACOBI, PETSCDTNODES_EQUISPACED, PETSCDTNODES_TANHSINH} PetscDTNodeType;
50d4afb720SToby Isaac 
51d4afb720SToby Isaac PETSC_EXTERN const char *const PetscDTNodeTypes[];
52d4afb720SToby Isaac 
5321454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *);
54c9638911SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *);
55bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt*);
56bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt);
57a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt*);
58a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt);
594f9ab2b4SJed Brown PETSC_EXTERN PetscErrorCode PetscQuadratureEqual(PetscQuadrature, PetscQuadrature, PetscBool*);
60a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt*, PetscInt*, PetscInt*, const PetscReal *[], const PetscReal *[]);
61a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal [], const PetscReal []);
6221454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer);
6321454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *);
64a0845e3aSMatthew G. Knepley 
652df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTensorQuadratureCreate(PetscQuadrature, PetscQuadrature, PetscQuadrature *);
6689710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *);
6789710940SMatthew G. Knepley 
68907761f8SToby Isaac PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *);
69907761f8SToby Isaac 
7037045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*);
71fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiNorm(PetscReal,PetscReal,PetscInt,PetscReal *);
7294e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt,PetscReal,PetscReal,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*);
73fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEvalJet(PetscReal,PetscReal,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscReal[]);
74fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPKDEvalJet(PetscInt,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscReal[]);
75d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedSize(PetscInt,PetscInt,PetscInt,PetscInt*);
76d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscInt,PetscReal[]);
7737045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt,PetscReal,PetscReal,PetscReal*,PetscReal*);
7894e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*,PetscReal*);
7994e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*,PetscReal*);
80916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt,PetscGaussLobattoLegendreCreateType,PetscReal*,PetscReal*);
81194825f6SJed Brown PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt,PetscInt,const PetscReal*,PetscInt,const PetscReal*,PetscReal*);
82a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
83e6a796c3SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
8437045ce4SJed Brown 
85b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
86d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
87d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
88b3c0f97bSTom Klotz 
89916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *);
90916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
91916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
92916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
93916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
94916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
95916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
96916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
97916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
98916e780bShannah_mairs 
991a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1001a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1011a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
1021a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
1031a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
1041a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1051a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
106dda711d0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
1071a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
1081a989b97SToby Isaac 
109d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt,PetscInt,const PetscInt[],PetscInt*);
110d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt,PetscInt,PetscInt,PetscInt[]);
111fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGradedOrderToIndex(PetscInt,const PetscInt[],PetscInt*);
112fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToGradedOrder(PetscInt,PetscInt,PetscInt[]);
113d4afb720SToby Isaac 
114fad4db65SToby Isaac #if defined(PETSC_USE_64BIT_INDICES)
115fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 20
116fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX  61
117fad4db65SToby Isaac #else
118fad4db65SToby Isaac #define PETSC_FACTORIAL_MAX 12
119fad4db65SToby Isaac #define PETSC_BINOMIAL_MAX  29
120fad4db65SToby Isaac #endif
121fad4db65SToby Isaac 
122fad4db65SToby Isaac /*MC
123fad4db65SToby Isaac    PetscDTFactorial - Approximate n! as a real number
124fad4db65SToby Isaac 
1254165533cSJose E. Roman    Input Parameter:
126fad4db65SToby Isaac .  n - a non-negative integer
127fad4db65SToby Isaac 
1284165533cSJose E. Roman    Output Parameter:
129fad4db65SToby Isaac .  factorial - n!
130fad4db65SToby Isaac 
131fad4db65SToby Isaac    Level: beginner
132fad4db65SToby Isaac M*/
1339fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial)
134fad4db65SToby Isaac {
135fad4db65SToby Isaac   PetscReal f = 1.0;
136fad4db65SToby Isaac 
137fad4db65SToby Isaac   PetscFunctionBegin;
138e2ab39ccSLisandro Dalcin   *factorial = -1.0;
13963a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %" PetscInt_FMT, n);
1405f80ce2aSJacob Faibussowitsch   for (PetscInt i = 1; i < n+1; ++i) f *= (PetscReal)i;
141fad4db65SToby Isaac   *factorial = f;
142fad4db65SToby Isaac   PetscFunctionReturn(0);
143fad4db65SToby Isaac }
144fad4db65SToby Isaac 
145fad4db65SToby Isaac /*MC
146fad4db65SToby Isaac    PetscDTFactorialInt - Compute n! as an integer
147fad4db65SToby Isaac 
1484165533cSJose E. Roman    Input Parameter:
149fad4db65SToby Isaac .  n - a non-negative integer
150fad4db65SToby Isaac 
1514165533cSJose E. Roman    Output Parameter:
152fad4db65SToby Isaac .  factorial - n!
153fad4db65SToby Isaac 
154fad4db65SToby Isaac    Level: beginner
155fad4db65SToby Isaac 
156fad4db65SToby 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.
157fad4db65SToby Isaac M*/
1589fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial)
159fad4db65SToby Isaac {
160fad4db65SToby Isaac   PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600};
161fad4db65SToby Isaac 
16228222859SToby Isaac   PetscFunctionBegin;
16328222859SToby Isaac   *factorial = -1;
16463a3b9bcSJacob 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);
165fad4db65SToby Isaac   if (n <= 12) {
166fad4db65SToby Isaac     *factorial = facLookup[n];
167fad4db65SToby Isaac   } else {
168fad4db65SToby Isaac     PetscInt f = facLookup[12];
169fad4db65SToby Isaac     PetscInt i;
170fad4db65SToby Isaac 
171fad4db65SToby Isaac     for (i = 13; i < n+1; ++i) f *= i;
172fad4db65SToby Isaac     *factorial = f;
173fad4db65SToby Isaac   }
174fad4db65SToby Isaac   PetscFunctionReturn(0);
175fad4db65SToby Isaac }
176fad4db65SToby Isaac 
177fad4db65SToby Isaac /*MC
178fad4db65SToby Isaac    PetscDTBinomial - Approximate the binomial coefficient "n choose k"
179fad4db65SToby Isaac 
1804165533cSJose E. Roman    Input Parameters:
181fad4db65SToby Isaac +  n - a non-negative integer
182fad4db65SToby Isaac -  k - an integer between 0 and n, inclusive
183fad4db65SToby Isaac 
1844165533cSJose E. Roman    Output Parameter:
185fad4db65SToby Isaac .  binomial - approximation of the binomial coefficient n choose k
186fad4db65SToby Isaac 
187fad4db65SToby Isaac    Level: beginner
188fad4db65SToby Isaac M*/
1899fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
1901a989b97SToby Isaac {
1911a989b97SToby Isaac   PetscFunctionBeginHot;
192e2ab39ccSLisandro Dalcin   *binomial = -1.0;
19363a3b9bcSJacob 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);
1941a989b97SToby Isaac   if (n <= 3) {
1951a989b97SToby Isaac     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};
1961a989b97SToby Isaac 
197e2ab39ccSLisandro Dalcin     *binomial = (PetscReal)binomLookup[n][k];
1981a989b97SToby Isaac   } else {
199e2ab39ccSLisandro Dalcin     PetscReal binom = 1.0;
2001a989b97SToby Isaac 
2011a989b97SToby Isaac     k = PetscMin(k, n - k);
2025f80ce2aSJacob Faibussowitsch     for (PetscInt 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;
22863a3b9bcSJacob 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);
22963a3b9bcSJacob 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);
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 
237fad4db65SToby Isaac     k = PetscMin(k, n - k);
2385f80ce2aSJacob Faibussowitsch     for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
23928222859SToby Isaac     bin = binom;
240fad4db65SToby Isaac   }
24128222859SToby Isaac   *binomial = bin;
242fad4db65SToby Isaac   PetscFunctionReturn(0);
243fad4db65SToby Isaac }
244fad4db65SToby Isaac 
245fad4db65SToby Isaac /*MC
246fad4db65SToby Isaac    PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps.
247fad4db65SToby Isaac 
248fad4db65SToby Isaac    A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation,
249fad4db65SToby Isaac    by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in
25028222859SToby Isaac    some position j >= i.  This swap is encoded as the difference (j - i).  The difference d_i at step i is less than
25128222859SToby Isaac    (n - i).  This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number
252fad4db65SToby Isaac    (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}.
253fad4db65SToby Isaac 
2544165533cSJose E. Roman    Input Parameters:
255fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
2568cd1e013SToby Isaac -  k - an integer in [0, n!)
257fad4db65SToby Isaac 
2584165533cSJose E. Roman    Output Parameters:
259fad4db65SToby Isaac +  perm - the permuted list of the integers [0, ..., n-1]
2608cd1e013SToby Isaac -  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.
261fad4db65SToby Isaac 
262fad4db65SToby 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.
263fad4db65SToby Isaac 
264fad4db65SToby Isaac    Level: beginner
265fad4db65SToby Isaac M*/
2669fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd)
2671a989b97SToby Isaac {
2681a989b97SToby Isaac   PetscInt  odd = 0;
2691a989b97SToby Isaac   PetscInt  i;
270fad4db65SToby Isaac   PetscInt  work[PETSC_FACTORIAL_MAX];
271fad4db65SToby Isaac   PetscInt *w;
2721a989b97SToby Isaac 
27328222859SToby Isaac   PetscFunctionBegin;
27428222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
27563a3b9bcSJacob 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);
276fad4db65SToby Isaac   w = &work[n - 2];
2771a989b97SToby Isaac   for (i = 2; i <= n; i++) {
2781a989b97SToby Isaac     *(w--) = k % i;
2791a989b97SToby Isaac     k /= i;
2801a989b97SToby Isaac   }
2811a989b97SToby Isaac   for (i = 0; i < n; i++) perm[i] = i;
2821a989b97SToby Isaac   for (i = 0; i < n - 1; i++) {
2831a989b97SToby Isaac     PetscInt s = work[i];
2841a989b97SToby Isaac     PetscInt swap = perm[i];
2851a989b97SToby Isaac 
2861a989b97SToby Isaac     perm[i] = perm[i + s];
2871a989b97SToby Isaac     perm[i + s] = swap;
2881a989b97SToby Isaac     odd ^= (!!s);
2891a989b97SToby Isaac   }
2901a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
2911a989b97SToby Isaac   PetscFunctionReturn(0);
2921a989b97SToby Isaac }
2931a989b97SToby Isaac 
294fad4db65SToby Isaac /*MC
2958cd1e013SToby Isaac    PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!).  This inverts PetscDTEnumPerm.
2968cd1e013SToby Isaac 
2974165533cSJose E. Roman    Input Parameters:
2988cd1e013SToby Isaac +  n - a non-negative integer (see note about limits below)
2998cd1e013SToby Isaac -  perm - the permuted list of the integers [0, ..., n-1]
3008cd1e013SToby Isaac 
3014165533cSJose E. Roman    Output Parameters:
3028cd1e013SToby Isaac +  k - an integer in [0, n!)
303f0fc11ceSJed Brown -  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.
3048cd1e013SToby Isaac 
3058cd1e013SToby 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.
3068cd1e013SToby Isaac 
3078cd1e013SToby Isaac    Level: beginner
3088cd1e013SToby Isaac M*/
3099fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd)
3108cd1e013SToby Isaac {
3118cd1e013SToby Isaac   PetscInt  odd = 0;
3128cd1e013SToby Isaac   PetscInt  i, idx;
3138cd1e013SToby Isaac   PetscInt  work[PETSC_FACTORIAL_MAX];
3148cd1e013SToby Isaac   PetscInt  iwork[PETSC_FACTORIAL_MAX];
3158cd1e013SToby Isaac 
3168cd1e013SToby Isaac   PetscFunctionBeginHot;
31728222859SToby Isaac   *k = -1;
31828222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
31963a3b9bcSJacob 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);
3208cd1e013SToby Isaac   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
3218cd1e013SToby Isaac   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
3228cd1e013SToby Isaac   for (idx = 0, i = 0; i < n - 1; i++) {
3238cd1e013SToby Isaac     PetscInt j = perm[i];
3248cd1e013SToby Isaac     PetscInt icur = work[i];
3258cd1e013SToby Isaac     PetscInt jloc = iwork[j];
3268cd1e013SToby Isaac     PetscInt diff = jloc - i;
3278cd1e013SToby Isaac 
3288cd1e013SToby Isaac     idx = idx * (n - i) + diff;
3298cd1e013SToby Isaac     /* swap (i, jloc) */
3308cd1e013SToby Isaac     work[i] = j;
3318cd1e013SToby Isaac     work[jloc] = icur;
3328cd1e013SToby Isaac     iwork[j] = i;
3338cd1e013SToby Isaac     iwork[icur] = jloc;
3348cd1e013SToby Isaac     odd ^= (!!diff);
3358cd1e013SToby Isaac   }
3368cd1e013SToby Isaac   *k = idx;
3378cd1e013SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
3388cd1e013SToby Isaac   PetscFunctionReturn(0);
3398cd1e013SToby Isaac }
3408cd1e013SToby Isaac 
3418cd1e013SToby Isaac /*MC
342fad4db65SToby Isaac    PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k).
343fad4db65SToby Isaac    The encoding is in lexicographic order.
344fad4db65SToby Isaac 
3454165533cSJose E. Roman    Input Parameters:
346fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
347fad4db65SToby Isaac .  k - an integer in [0, n]
348fad4db65SToby Isaac -  j - an index in [0, n choose k)
349fad4db65SToby Isaac 
3504165533cSJose E. Roman    Output Parameter:
351fad4db65SToby Isaac .  subset - the jth subset of size k of the integers [0, ..., n - 1]
352fad4db65SToby Isaac 
353fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
354fad4db65SToby Isaac 
355fad4db65SToby Isaac    Level: beginner
356fad4db65SToby Isaac 
357db781477SPatrick Sanan .seealso: `PetscDTSubsetIndex()`
358fad4db65SToby Isaac M*/
3599fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
3601a989b97SToby Isaac {
3615f80ce2aSJacob Faibussowitsch   PetscInt Nk;
3621a989b97SToby Isaac 
3631a989b97SToby Isaac   PetscFunctionBeginHot;
3649566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
3655f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
3661a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
3671a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
3681a989b97SToby Isaac 
3691a989b97SToby Isaac     if (j < Nminuskminus) {
3701a989b97SToby Isaac       subset[l++] = i;
3711a989b97SToby Isaac       Nk = Nminuskminus;
3721a989b97SToby Isaac     } else {
3731a989b97SToby Isaac       j -= Nminuskminus;
3741a989b97SToby Isaac       Nk = Nminusk;
3751a989b97SToby Isaac     }
3761a989b97SToby Isaac   }
3771a989b97SToby Isaac   PetscFunctionReturn(0);
3781a989b97SToby Isaac }
3791a989b97SToby Isaac 
380fad4db65SToby Isaac /*MC
381fad4db65SToby 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.
382fad4db65SToby Isaac 
3834165533cSJose E. Roman    Input Parameters:
384fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
385fad4db65SToby Isaac .  k - an integer in [0, n]
386fad4db65SToby Isaac -  subset - an ordered subset of the integers [0, ..., n - 1]
387fad4db65SToby Isaac 
3884165533cSJose E. Roman    Output Parameter:
389fad4db65SToby Isaac .  index - the rank of the subset in lexicographic order
390fad4db65SToby Isaac 
391fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
392fad4db65SToby Isaac 
393fad4db65SToby Isaac    Level: beginner
394fad4db65SToby Isaac 
395db781477SPatrick Sanan .seealso: `PetscDTEnumSubset()`
396fad4db65SToby Isaac M*/
3979fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
3981a989b97SToby Isaac {
3995f80ce2aSJacob Faibussowitsch   PetscInt j = 0, Nk;
4001a989b97SToby Isaac 
40128222859SToby Isaac   PetscFunctionBegin;
40228222859SToby Isaac   *index = -1;
4039566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
4045f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
4051a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
4061a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
4071a989b97SToby Isaac 
4081a989b97SToby Isaac     if (subset[l] == i) {
4091a989b97SToby Isaac       l++;
4101a989b97SToby Isaac       Nk = Nminuskminus;
4111a989b97SToby Isaac     } else {
4121a989b97SToby Isaac       j += Nminuskminus;
4131a989b97SToby Isaac       Nk = Nminusk;
4141a989b97SToby Isaac     }
4151a989b97SToby Isaac   }
4161a989b97SToby Isaac   *index = j;
4171a989b97SToby Isaac   PetscFunctionReturn(0);
4181a989b97SToby Isaac }
4191a989b97SToby Isaac 
420fad4db65SToby Isaac /*MC
42128222859SToby 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.
422fad4db65SToby Isaac 
4234165533cSJose E. Roman    Input Parameters:
424fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
425fad4db65SToby Isaac .  k - an integer in [0, n]
426fad4db65SToby Isaac -  j - an index in [0, n choose k)
427fad4db65SToby Isaac 
4284165533cSJose E. Roman    Output Parameters:
429fad4db65SToby Isaac +  perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set.
43028222859SToby Isaac -  isOdd - if not NULL, return whether perm is an even or odd permutation.
431fad4db65SToby Isaac 
432fad4db65SToby Isaac    Note: this is limited by arguments such that n choose k can be represented by PetscInt
433fad4db65SToby Isaac 
434fad4db65SToby Isaac    Level: beginner
435fad4db65SToby Isaac 
436db781477SPatrick Sanan .seealso: `PetscDTEnumSubset()`, `PetscDTSubsetIndex()`
437fad4db65SToby Isaac M*/
4389fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd)
4391a989b97SToby Isaac {
4405f80ce2aSJacob Faibussowitsch   PetscInt i, l, m, Nk, odd = 0;
4415f80ce2aSJacob Faibussowitsch   PetscInt *subcomp = perm+k;
4421a989b97SToby Isaac 
44328222859SToby Isaac   PetscFunctionBegin;
44428222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
4459566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
4461a989b97SToby Isaac   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
4471a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
4481a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
4491a989b97SToby Isaac 
4501a989b97SToby Isaac     if (j < Nminuskminus) {
451fad4db65SToby Isaac       perm[l++] = i;
4521a989b97SToby Isaac       Nk = Nminuskminus;
4531a989b97SToby Isaac     } else {
4541a989b97SToby Isaac       subcomp[m++] = i;
4551a989b97SToby Isaac       j -= Nminuskminus;
4561a989b97SToby Isaac       odd ^= ((k - l) & 1);
4571a989b97SToby Isaac       Nk = Nminusk;
4581a989b97SToby Isaac     }
4591a989b97SToby Isaac   }
4605f80ce2aSJacob Faibussowitsch   for (; i < n; i++) subcomp[m++] = i;
4611a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
4621a989b97SToby Isaac   PetscFunctionReturn(0);
4631a989b97SToby Isaac }
4641a989b97SToby Isaac 
465ef0bb6c7SMatthew G. Knepley struct _p_PetscTabulation {
466a5b23f4aSJose E. Roman   PetscInt    K;    /* Indicates a k-jet, namely tabulated derivatives up to order k */
46719815104SMartin Diehl   PetscInt    Nr;   /* The number of tabulation replicas (often 1) */
468ef0bb6c7SMatthew G. Knepley   PetscInt    Np;   /* The number of tabulation points in a replica */
469ef0bb6c7SMatthew G. Knepley   PetscInt    Nb;   /* The number of functions tabulated */
470ef0bb6c7SMatthew G. Knepley   PetscInt    Nc;   /* The number of function components */
471ef0bb6c7SMatthew G. Knepley   PetscInt    cdim; /* The coordinate dimension */
472ef0bb6c7SMatthew G. Knepley   PetscReal **T;    /* The tabulation T[K] of functions and their derivatives
473ef0bb6c7SMatthew G. Knepley                        T[0] = B[Nr*Np][Nb][Nc]:             The basis function values at quadrature points
474ef0bb6c7SMatthew G. Knepley                        T[1] = D[Nr*Np][Nb][Nc][cdim]:       The basis function derivatives at quadrature points
475ef0bb6c7SMatthew G. Knepley                        T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */
476ef0bb6c7SMatthew G. Knepley };
477ef0bb6c7SMatthew G. Knepley typedef struct _p_PetscTabulation *PetscTabulation;
478ef0bb6c7SMatthew G. Knepley 
479d6685f55SMatthew G. Knepley typedef PetscErrorCode (*PetscProbFunc)(const PetscReal[], const PetscReal[], PetscReal[]);
480d6685f55SMatthew G. Knepley 
481d6685f55SMatthew G. Knepley typedef enum {DTPROB_DENSITY_CONSTANT, DTPROB_DENSITY_GAUSSIAN, DTPROB_DENSITY_MAXWELL_BOLTZMANN, DTPROB_NUM_DENSITY} DTProbDensityType;
482d6685f55SMatthew G. Knepley PETSC_EXTERN const char * const DTProbDensityTypes[];
483d6685f55SMatthew G. Knepley 
484d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
485d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
486d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
487d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
488d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
489d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
490d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
491d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
492d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
493d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
494d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
495ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]);
496ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]);
497d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
498d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
499d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
500ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
501ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
502ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
503ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
504ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscCDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
505ea1b28ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant3D(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