xref: /petsc/include/petscdt.h (revision daa252a5865999f1c6bc8348e4bad0b69e52d518)
137045ce4SJed Brown /*
237045ce4SJed Brown   Common tools for constructing discretizations
337045ce4SJed Brown */
4a4963045SJacob Faibussowitsch #pragma once
537045ce4SJed Brown 
637045ce4SJed Brown #include <petscsys.h>
74366bac7SMatthew G. Knepley #include <petscdmtypes.h>
84366bac7SMatthew G. Knepley #include <petscistypes.h>
937045ce4SJed Brown 
10ce78bad3SBarry Smith /* MANSEC = DM */
11ac09b921SBarry Smith /* SUBMANSEC = DT */
12ac09b921SBarry Smith 
132cd22861SMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCQUADRATURE_CLASSID;
142cd22861SMatthew G. Knepley 
1521454ff5SMatthew G. Knepley /*S
1616a05f60SBarry Smith   PetscQuadrature - Quadrature rule for numerical integration.
1721454ff5SMatthew G. Knepley 
18329bbf4eSMatthew G. Knepley   Level: beginner
1921454ff5SMatthew G. Knepley 
20db781477SPatrick Sanan .seealso: `PetscQuadratureCreate()`, `PetscQuadratureDestroy()`
2121454ff5SMatthew G. Knepley S*/
2221454ff5SMatthew G. Knepley typedef struct _p_PetscQuadrature *PetscQuadrature;
2321454ff5SMatthew G. Knepley 
248272889dSSatish Balay /*E
25916e780bShannah_mairs   PetscGaussLobattoLegendreCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights
268272889dSSatish Balay 
2716a05f60SBarry Smith   Values:
2816a05f60SBarry Smith +  `PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA` - compute the nodes via linear algebra
2916a05f60SBarry Smith -  `PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON`         - compute the nodes by solving a nonlinear equation with Newton's method
3016a05f60SBarry Smith 
318272889dSSatish Balay   Level: intermediate
328272889dSSatish Balay 
3316a05f60SBarry Smith .seealso: `PetscQuadrature`
348272889dSSatish Balay E*/
359371c9d4SSatish Balay typedef enum {
369371c9d4SSatish Balay   PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,
379371c9d4SSatish Balay   PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON
389371c9d4SSatish Balay } PetscGaussLobattoLegendreCreateType;
398272889dSSatish Balay 
40d4afb720SToby Isaac /*E
41d4afb720SToby Isaac   PetscDTNodeType - A description of strategies for generating nodes (both
42d4afb720SToby Isaac   quadrature nodes and nodes for Lagrange polynomials)
43d4afb720SToby Isaac 
4416a05f60SBarry Smith   Values:
4516a05f60SBarry Smith + `PETSCDTNODES_DEFAULT`     - Nodes chosen by PETSc
4616a05f60SBarry Smith . `PETSCDTNODES_GAUSSJACOBI` - Nodes at either Gauss-Jacobi or Gauss-Lobatto-Jacobi quadrature points
4716a05f60SBarry Smith . `PETSCDTNODES_EQUISPACED`  - Nodes equispaced either including the endpoints or excluding them
4816a05f60SBarry Smith - `PETSCDTNODES_TANHSINH`    - Nodes at Tanh-Sinh quadrature points
49d4afb720SToby Isaac 
5016a05f60SBarry Smith   Level: intermediate
51d4afb720SToby Isaac 
5287497f52SBarry Smith   Note:
5387497f52SBarry Smith   A `PetscDTNodeType` can be paired with a `PetscBool` to indicate whether
5487497f52SBarry Smith   the nodes include endpoints or not, and in the case of `PETSCDT_GAUSSJACOBI`
55d4afb720SToby Isaac   with exponents for the weight function.
56d4afb720SToby Isaac 
5716a05f60SBarry Smith .seealso: `PetscQuadrature`
58d4afb720SToby Isaac E*/
599371c9d4SSatish Balay typedef enum {
609371c9d4SSatish Balay   PETSCDTNODES_DEFAULT     = -1,
61ce78bad3SBarry Smith   PETSCDTNODES_GAUSSJACOBI = 0,
62ce78bad3SBarry Smith   PETSCDTNODES_EQUISPACED  = 1,
63ce78bad3SBarry Smith   PETSCDTNODES_TANHSINH    = 2
649371c9d4SSatish Balay } PetscDTNodeType;
65d4afb720SToby Isaac 
66d3c69ad0SToby Isaac PETSC_EXTERN const char *const *const PetscDTNodeTypes;
67d3c69ad0SToby Isaac 
68d3c69ad0SToby Isaac /*E
69d3c69ad0SToby Isaac   PetscDTSimplexQuadratureType - A description of classes of quadrature rules for simplices
70d3c69ad0SToby Isaac 
7116a05f60SBarry Smith   Values:
7216a05f60SBarry Smith +  `PETSCDTSIMPLEXQUAD_DEFAULT` - Quadrature rule chosen by PETSc
7316a05f60SBarry Smith .  `PETSCDTSIMPLEXQUAD_CONIC`   - Quadrature rules constructed as
74d3c69ad0SToby Isaac                                   conically-warped tensor products of 1D
75d3c69ad0SToby Isaac                                   Gauss-Jacobi quadrature rules.  These are
76d3c69ad0SToby Isaac                                   explicitly computable in any dimension for any
77d3c69ad0SToby Isaac                                   degree, and the tensor-product structure can be
78d3c69ad0SToby Isaac                                   exploited by sum-factorization methods, but
79d3c69ad0SToby Isaac                                   they are not efficient in terms of nodes per
80d3c69ad0SToby Isaac                                   polynomial degree.
8116a05f60SBarry Smith -  `PETSCDTSIMPLEXQUAD_MINSYM`  - Quadrature rules that are fully symmetric
82d3c69ad0SToby Isaac                                   (symmetries of the simplex preserve the nodes
83d3c69ad0SToby Isaac                                   and weights) with minimal (or near minimal)
84d3c69ad0SToby Isaac                                   number of nodes.  In dimensions higher than 1
85d3c69ad0SToby Isaac                                   these are not simple to compute, so lookup
86d3c69ad0SToby Isaac                                   tables are used.
87d3c69ad0SToby Isaac 
8816a05f60SBarry Smith   Level: intermediate
8916a05f60SBarry Smith 
9016a05f60SBarry Smith .seealso: `PetscQuadrature`, `PetscDTSimplexQuadrature()`
91d3c69ad0SToby Isaac E*/
929371c9d4SSatish Balay typedef enum {
939371c9d4SSatish Balay   PETSCDTSIMPLEXQUAD_DEFAULT = -1,
949371c9d4SSatish Balay   PETSCDTSIMPLEXQUAD_CONIC   = 0,
95ce78bad3SBarry Smith   PETSCDTSIMPLEXQUAD_MINSYM  = 1
969371c9d4SSatish Balay } PetscDTSimplexQuadratureType;
97d3c69ad0SToby Isaac 
98d3c69ad0SToby Isaac PETSC_EXTERN const char *const *const PetscDTSimplexQuadratureTypes;
99d4afb720SToby Isaac 
10021454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *);
101c9638911SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *);
1024366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetCellType(PetscQuadrature, DMPolytopeType *);
1034366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetCellType(PetscQuadrature, DMPolytopeType);
104bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt *);
105bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt);
106a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt *);
107a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt);
1084f9ab2b4SJed Brown PETSC_EXTERN PetscErrorCode PetscQuadratureEqual(PetscQuadrature, PetscQuadrature, PetscBool *);
109a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt *, PetscInt *, PetscInt *, const PetscReal *[], const PetscReal *[]);
110a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal[], const PetscReal[]);
11121454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer);
11221454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *);
113a0845e3aSMatthew G. Knepley 
1142df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTensorQuadratureCreate(PetscQuadrature, PetscQuadrature, PetscQuadrature *);
11589710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *);
1164366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureComputePermutations(PetscQuadrature, PetscInt *, IS *[]);
11789710940SMatthew G. Knepley 
118907761f8SToby Isaac PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *);
119907761f8SToby Isaac 
12037045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *);
121fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiNorm(PetscReal, PetscReal, PetscInt, PetscReal *);
12294e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt, PetscReal, PetscReal, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *);
123fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTJacobiEvalJet(PetscReal, PetscReal, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]);
124fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPKDEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]);
125d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedSize(PetscInt, PetscInt, PetscInt, PetscInt *);
126d8f25ad8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscInt, PetscReal[]);
12737045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt, PetscReal, PetscReal, PetscReal *, PetscReal *);
12894e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *);
12994e21283SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *);
130916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt, PetscGaussLobattoLegendreCreateType, PetscReal *, PetscReal *);
131194825f6SJed Brown PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
132a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
133e6a796c3SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
134d3c69ad0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTSimplexQuadrature(PetscInt, PetscInt, PetscDTSimplexQuadratureType, PetscQuadrature *);
1354366bac7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTCreateDefaultQuadrature(DMPolytopeType, PetscInt, PetscQuadrature *, PetscQuadrature *);
136*f2c64c88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTCreateQuadratureByCell(DMPolytopeType, PetscInt, PetscDTSimplexQuadratureType, PetscQuadrature *, PetscQuadrature *);
13737045ce4SJed Brown 
138b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
139d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
140d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
141b3c0f97bSTom Klotz 
142916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *);
143916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
144916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
145916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
146916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
147916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
148916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
149916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
150916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
151916e780bShannah_mairs 
1522dce792eSToby Isaac /*MC
1532dce792eSToby Isaac   PETSC_FORM_DEGREE_UNDEFINED - Indicates that a field does not have
1542dce792eSToby Isaac   a well-defined form degree in exterior calculus.
1552dce792eSToby Isaac 
1562dce792eSToby Isaac   Level: advanced
1572dce792eSToby Isaac 
1582dce792eSToby Isaac .seealso: `PetscDTAltV`, `PetscDualSpaceGetFormDegree()`
1592dce792eSToby Isaac M*/
1602dce792eSToby Isaac #define PETSC_FORM_DEGREE_UNDEFINED PETSC_INT_MIN
1612dce792eSToby Isaac 
1621a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1631a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1641a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
1651a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
1661a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
1671a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1681a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
169dda711d0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
1701a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
1711a989b97SToby Isaac 
172d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt, PetscInt, const PetscInt[], PetscInt *);
173d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt, PetscInt, PetscInt, PetscInt[]);
174fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGradedOrderToIndex(PetscInt, const PetscInt[], PetscInt *);
175fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToGradedOrder(PetscInt, PetscInt, PetscInt[]);
176d4afb720SToby Isaac 
177fad4db65SToby Isaac #if defined(PETSC_USE_64BIT_INDICES)
178fad4db65SToby Isaac   #define PETSC_FACTORIAL_MAX 20
179fad4db65SToby Isaac   #define PETSC_BINOMIAL_MAX  61
180fad4db65SToby Isaac #else
181fad4db65SToby Isaac   #define PETSC_FACTORIAL_MAX 12
182fad4db65SToby Isaac   #define PETSC_BINOMIAL_MAX  29
183fad4db65SToby Isaac #endif
184fad4db65SToby Isaac 
185fad4db65SToby Isaac /*MC
186fad4db65SToby Isaac    PetscDTFactorial - Approximate n! as a real number
187fad4db65SToby Isaac 
1884165533cSJose E. Roman    Input Parameter:
189fad4db65SToby Isaac .  n - a non-negative integer
190fad4db65SToby Isaac 
1914165533cSJose E. Roman    Output Parameter:
192fad4db65SToby Isaac .  factorial - n!
193fad4db65SToby Isaac 
194fad4db65SToby Isaac    Level: beginner
19516a05f60SBarry Smith 
19616a05f60SBarry Smith .seealso: `PetscDTFactorialInt()`, `PetscDTBinomialInt()`, `PetscDTBinomial()`
197fad4db65SToby Isaac M*/
PetscDTFactorial(PetscInt n,PetscReal * factorial)198d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial)
199d71ae5a4SJacob Faibussowitsch {
200fad4db65SToby Isaac   PetscReal f = 1.0;
201fad4db65SToby Isaac 
202fad4db65SToby Isaac   PetscFunctionBegin;
203e2ab39ccSLisandro Dalcin   *factorial = -1.0;
20463a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %" PetscInt_FMT, n);
2055f80ce2aSJacob Faibussowitsch   for (PetscInt i = 1; i < n + 1; ++i) f *= (PetscReal)i;
206fad4db65SToby Isaac   *factorial = f;
2073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
208fad4db65SToby Isaac }
209fad4db65SToby Isaac 
210fad4db65SToby Isaac /*MC
211fad4db65SToby Isaac    PetscDTFactorialInt - Compute n! as an integer
212fad4db65SToby Isaac 
2134165533cSJose E. Roman    Input Parameter:
214fad4db65SToby Isaac .  n - a non-negative integer
215fad4db65SToby Isaac 
2164165533cSJose E. Roman    Output Parameter:
217fad4db65SToby Isaac .  factorial - n!
218fad4db65SToby Isaac 
219fad4db65SToby Isaac    Level: beginner
220fad4db65SToby Isaac 
22116a05f60SBarry Smith    Note:
22295bd0b28SBarry Smith    This is limited to `n` such that n! can be represented by `PetscInt`, which is 12 if `PetscInt` is a signed 32-bit integer and 20 if `PetscInt` is a signed 64-bit integer.
22316a05f60SBarry Smith 
22416a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTBinomialInt()`, `PetscDTBinomial()`
225fad4db65SToby Isaac M*/
PetscDTFactorialInt(PetscInt n,PetscInt * factorial)226d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial)
227d71ae5a4SJacob Faibussowitsch {
228fad4db65SToby Isaac   PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600};
229fad4db65SToby Isaac 
23028222859SToby Isaac   PetscFunctionBegin;
23128222859SToby Isaac   *factorial = -1;
23263a3b9bcSJacob 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);
233fad4db65SToby Isaac   if (n <= 12) {
234fad4db65SToby Isaac     *factorial = facLookup[n];
235fad4db65SToby Isaac   } else {
236fad4db65SToby Isaac     PetscInt f = facLookup[12];
237fad4db65SToby Isaac     PetscInt i;
238fad4db65SToby Isaac 
239fad4db65SToby Isaac     for (i = 13; i < n + 1; ++i) f *= i;
240fad4db65SToby Isaac     *factorial = f;
241fad4db65SToby Isaac   }
2423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
243fad4db65SToby Isaac }
244fad4db65SToby Isaac 
245fad4db65SToby Isaac /*MC
24695bd0b28SBarry Smith    PetscDTBinomial - Approximate the binomial coefficient `n` choose `k`
247fad4db65SToby Isaac 
2484165533cSJose E. Roman    Input Parameters:
249fad4db65SToby Isaac +  n - a non-negative integer
25095bd0b28SBarry Smith -  k - an integer between 0 and `n`, inclusive
251fad4db65SToby Isaac 
2524165533cSJose E. Roman    Output Parameter:
25395bd0b28SBarry Smith .  binomial - approximation of the binomial coefficient `n` choose `k`
254fad4db65SToby Isaac 
255fad4db65SToby Isaac    Level: beginner
25616a05f60SBarry Smith 
25716a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`
258fad4db65SToby Isaac M*/
PetscDTBinomial(PetscInt n,PetscInt k,PetscReal * binomial)259d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
260d71ae5a4SJacob Faibussowitsch {
2611a989b97SToby Isaac   PetscFunctionBeginHot;
262e2ab39ccSLisandro Dalcin   *binomial = -1.0;
26363a3b9bcSJacob 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);
2641a989b97SToby Isaac   if (n <= 3) {
2659371c9d4SSatish Balay     PetscInt binomLookup[4][4] = {
2669371c9d4SSatish Balay       {1, 0, 0, 0},
2679371c9d4SSatish Balay       {1, 1, 0, 0},
2689371c9d4SSatish Balay       {1, 2, 1, 0},
2699371c9d4SSatish Balay       {1, 3, 3, 1}
2709371c9d4SSatish Balay     };
2711a989b97SToby Isaac 
272e2ab39ccSLisandro Dalcin     *binomial = (PetscReal)binomLookup[n][k];
2731a989b97SToby Isaac   } else {
274e2ab39ccSLisandro Dalcin     PetscReal binom = 1.0;
2751a989b97SToby Isaac 
2761a989b97SToby Isaac     k = PetscMin(k, n - k);
2775f80ce2aSJacob Faibussowitsch     for (PetscInt i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1);
2781a989b97SToby Isaac     *binomial = binom;
2791a989b97SToby Isaac   }
2803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2811a989b97SToby Isaac }
2821a989b97SToby Isaac 
283fad4db65SToby Isaac /*MC
28495bd0b28SBarry Smith    PetscDTBinomialInt - Compute the binomial coefficient `n` choose `k`
285fad4db65SToby Isaac 
28697bb3fdcSJose E. Roman    Input Parameters:
287fad4db65SToby Isaac +  n - a non-negative integer
28895bd0b28SBarry Smith -  k - an integer between 0 and `n`, inclusive
289fad4db65SToby Isaac 
29097bb3fdcSJose E. Roman    Output Parameter:
29195bd0b28SBarry Smith .  binomial - the binomial coefficient `n` choose `k`
292fad4db65SToby Isaac 
293fad4db65SToby Isaac    Level: beginner
29416a05f60SBarry Smith 
29516a05f60SBarry Smith    Note:
29616a05f60SBarry Smith    This is limited by integers that can be represented by `PetscInt`.
29716a05f60SBarry Smith 
29816a05f60SBarry Smith    Use `PetscDTBinomial()` for real number approximations of larger values
29916a05f60SBarry Smith 
30016a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTEnumPerm()`
301fad4db65SToby Isaac M*/
PetscDTBinomialInt(PetscInt n,PetscInt k,PetscInt * binomial)302d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial)
303d71ae5a4SJacob Faibussowitsch {
30428222859SToby Isaac   PetscInt bin;
30528222859SToby Isaac 
30628222859SToby Isaac   PetscFunctionBegin;
30728222859SToby Isaac   *binomial = -1;
30863a3b9bcSJacob 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);
30963a3b9bcSJacob 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);
310fad4db65SToby Isaac   if (n <= 3) {
3119371c9d4SSatish Balay     PetscInt binomLookup[4][4] = {
3129371c9d4SSatish Balay       {1, 0, 0, 0},
3139371c9d4SSatish Balay       {1, 1, 0, 0},
3149371c9d4SSatish Balay       {1, 2, 1, 0},
3159371c9d4SSatish Balay       {1, 3, 3, 1}
3169371c9d4SSatish Balay     };
317fad4db65SToby Isaac 
31828222859SToby Isaac     bin = binomLookup[n][k];
319fad4db65SToby Isaac   } else {
320fad4db65SToby Isaac     PetscInt binom = 1;
321fad4db65SToby Isaac 
322fad4db65SToby Isaac     k = PetscMin(k, n - k);
3235f80ce2aSJacob Faibussowitsch     for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
32428222859SToby Isaac     bin = binom;
325fad4db65SToby Isaac   }
32628222859SToby Isaac   *binomial = bin;
3273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
328fad4db65SToby Isaac }
329fad4db65SToby Isaac 
330ce78bad3SBarry Smith /* the following inline routines should be not be inline routines and then Fortran binding can be built automatically */
331ce78bad3SBarry Smith #define PeOp
332ce78bad3SBarry Smith 
333fad4db65SToby Isaac /*MC
33416a05f60SBarry Smith    PetscDTEnumPerm - Get a permutation of `n` integers from its encoding into the integers [0, n!) as a sequence of swaps.
335fad4db65SToby Isaac 
3364165533cSJose E. Roman    Input Parameters:
337fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
3388cd1e013SToby Isaac -  k - an integer in [0, n!)
339fad4db65SToby Isaac 
3404165533cSJose E. Roman    Output Parameters:
341fad4db65SToby Isaac +  perm  - the permuted list of the integers [0, ..., n-1]
3422fe279fdSBarry Smith -  isOdd - if not `NULL`, returns whether the permutation used an even or odd number of swaps.
343fad4db65SToby Isaac 
34416a05f60SBarry Smith    Level: intermediate
345fad4db65SToby Isaac 
34616a05f60SBarry Smith    Notes:
34716a05f60SBarry Smith    A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation,
34816a05f60SBarry Smith    by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in
34916a05f60SBarry Smith    some position j >= i.  This swap is encoded as the difference (j - i).  The difference d_i at step i is less than
35016a05f60SBarry Smith    (n - i).  This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number
35116a05f60SBarry Smith    (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}.
35216a05f60SBarry Smith 
35316a05f60SBarry Smith    Limited to `n` such that `n`! can be represented by `PetscInt`, which is 12 if `PetscInt` is a signed 32-bit integer and 20 if `PetscInt` is a signed 64-bit integer.
35416a05f60SBarry Smith 
35516a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTPermIndex()`
356fad4db65SToby Isaac M*/
PetscDTEnumPerm(PetscInt n,PetscInt k,PetscInt * perm,PeOp PetscBool * isOdd)357ce78bad3SBarry Smith static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PeOp PetscBool *isOdd)
358d71ae5a4SJacob Faibussowitsch {
3591a989b97SToby Isaac   PetscInt  odd = 0;
3601a989b97SToby Isaac   PetscInt  i;
361fad4db65SToby Isaac   PetscInt  work[PETSC_FACTORIAL_MAX];
362fad4db65SToby Isaac   PetscInt *w;
3631a989b97SToby Isaac 
36428222859SToby Isaac   PetscFunctionBegin;
36528222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
36663a3b9bcSJacob 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);
367810441c8SPierre Jolivet   if (n >= 2) {
368fad4db65SToby Isaac     w = &work[n - 2];
3691a989b97SToby Isaac     for (i = 2; i <= n; i++) {
3701a989b97SToby Isaac       *(w--) = k % i;
3711a989b97SToby Isaac       k /= i;
3721a989b97SToby Isaac     }
373810441c8SPierre Jolivet   }
3741a989b97SToby Isaac   for (i = 0; i < n; i++) perm[i] = i;
3751a989b97SToby Isaac   for (i = 0; i < n - 1; i++) {
3761a989b97SToby Isaac     PetscInt s    = work[i];
3771a989b97SToby Isaac     PetscInt swap = perm[i];
3781a989b97SToby Isaac 
3791a989b97SToby Isaac     perm[i]     = perm[i + s];
3801a989b97SToby Isaac     perm[i + s] = swap;
3811a989b97SToby Isaac     odd ^= (!!s);
3821a989b97SToby Isaac   }
3831a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
3843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3851a989b97SToby Isaac }
3861a989b97SToby Isaac 
387fad4db65SToby Isaac /*MC
38816a05f60SBarry Smith    PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!).  This inverts `PetscDTEnumPerm()`.
3898cd1e013SToby Isaac 
3904165533cSJose E. Roman    Input Parameters:
3918cd1e013SToby Isaac +  n    - a non-negative integer (see note about limits below)
3928cd1e013SToby Isaac -  perm - the permuted list of the integers [0, ..., n-1]
3938cd1e013SToby Isaac 
3944165533cSJose E. Roman    Output Parameters:
3958cd1e013SToby Isaac +  k     - an integer in [0, n!)
3962fe279fdSBarry Smith -  isOdd - if not `NULL`, returns whether the permutation used an even or odd number of swaps.
3978cd1e013SToby Isaac 
3988cd1e013SToby Isaac    Level: beginner
39916a05f60SBarry Smith 
40016a05f60SBarry Smith    Note:
40116a05f60SBarry Smith    Limited to `n` such that `n`! can be represented by `PetscInt`, which is 12 if `PetscInt` is a signed 32-bit integer and 20 if `PetscInt` is a signed 64-bit integer.
40216a05f60SBarry Smith 
40316a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`
4048cd1e013SToby Isaac M*/
PetscDTPermIndex(PetscInt n,const PetscInt * perm,PetscInt * k,PeOp PetscBool * isOdd)405ce78bad3SBarry Smith static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PeOp PetscBool *isOdd)
406d71ae5a4SJacob Faibussowitsch {
4078cd1e013SToby Isaac   PetscInt odd = 0;
4088cd1e013SToby Isaac   PetscInt i, idx;
4098cd1e013SToby Isaac   PetscInt work[PETSC_FACTORIAL_MAX];
4108cd1e013SToby Isaac   PetscInt iwork[PETSC_FACTORIAL_MAX];
4118cd1e013SToby Isaac 
4128cd1e013SToby Isaac   PetscFunctionBeginHot;
41328222859SToby Isaac   *k = -1;
41428222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
41563a3b9bcSJacob 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);
4168cd1e013SToby Isaac   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
4178cd1e013SToby Isaac   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
4188cd1e013SToby Isaac   for (idx = 0, i = 0; i < n - 1; i++) {
4198cd1e013SToby Isaac     PetscInt j    = perm[i];
4208cd1e013SToby Isaac     PetscInt icur = work[i];
4218cd1e013SToby Isaac     PetscInt jloc = iwork[j];
4228cd1e013SToby Isaac     PetscInt diff = jloc - i;
4238cd1e013SToby Isaac 
4248cd1e013SToby Isaac     idx = idx * (n - i) + diff;
4258cd1e013SToby Isaac     /* swap (i, jloc) */
4268cd1e013SToby Isaac     work[i]     = j;
4278cd1e013SToby Isaac     work[jloc]  = icur;
4288cd1e013SToby Isaac     iwork[j]    = i;
4298cd1e013SToby Isaac     iwork[icur] = jloc;
4308cd1e013SToby Isaac     odd ^= (!!diff);
4318cd1e013SToby Isaac   }
4328cd1e013SToby Isaac   *k = idx;
4338cd1e013SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
4343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4358cd1e013SToby Isaac }
4368cd1e013SToby Isaac 
4378cd1e013SToby Isaac /*MC
438fad4db65SToby Isaac    PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k).
439fad4db65SToby Isaac    The encoding is in lexicographic order.
440fad4db65SToby Isaac 
4414165533cSJose E. Roman    Input Parameters:
442fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
443fad4db65SToby Isaac .  k - an integer in [0, n]
444fad4db65SToby Isaac -  j - an index in [0, n choose k)
445fad4db65SToby Isaac 
4464165533cSJose E. Roman    Output Parameter:
447fad4db65SToby Isaac .  subset - the jth subset of size k of the integers [0, ..., n - 1]
448fad4db65SToby Isaac 
449fad4db65SToby Isaac    Level: beginner
450fad4db65SToby Isaac 
45116a05f60SBarry Smith    Note:
45216a05f60SBarry Smith    Limited by arguments such that `n` choose `k` can be represented by `PetscInt`
45316a05f60SBarry Smith 
45416a05f60SBarry Smith .seealso: `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()`
455fad4db65SToby Isaac M*/
PetscDTEnumSubset(PetscInt n,PetscInt k,PetscInt j,PetscInt * subset)456d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
457d71ae5a4SJacob Faibussowitsch {
4585f80ce2aSJacob Faibussowitsch   PetscInt Nk;
4591a989b97SToby Isaac 
4601a989b97SToby Isaac   PetscFunctionBeginHot;
4619566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
4625f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
4631a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
4641a989b97SToby Isaac     PetscInt Nminusk      = Nk - Nminuskminus;
4651a989b97SToby Isaac 
4661a989b97SToby Isaac     if (j < Nminuskminus) {
4671a989b97SToby Isaac       subset[l++] = i;
4681a989b97SToby Isaac       Nk          = Nminuskminus;
4691a989b97SToby Isaac     } else {
4701a989b97SToby Isaac       j -= Nminuskminus;
4711a989b97SToby Isaac       Nk = Nminusk;
4721a989b97SToby Isaac     }
4731a989b97SToby Isaac   }
4743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4751a989b97SToby Isaac }
4761a989b97SToby Isaac 
477fad4db65SToby Isaac /*MC
47887497f52SBarry Smith    PetscDTSubsetIndex - Convert an ordered subset of k integers from the set [0, ..., n - 1] to its encoding as an integers in [0, n choose k) in lexicographic order.
47987497f52SBarry Smith    This is the inverse of `PetscDTEnumSubset`.
480fad4db65SToby Isaac 
4814165533cSJose E. Roman    Input Parameters:
482fad4db65SToby Isaac +  n      - a non-negative integer (see note about limits below)
483fad4db65SToby Isaac .  k      - an integer in [0, n]
484fad4db65SToby Isaac -  subset - an ordered subset of the integers [0, ..., n - 1]
485fad4db65SToby Isaac 
4864165533cSJose E. Roman    Output Parameter:
487fad4db65SToby Isaac .  index - the rank of the subset in lexicographic order
488fad4db65SToby Isaac 
489fad4db65SToby Isaac    Level: beginner
490fad4db65SToby Isaac 
49116a05f60SBarry Smith    Note:
49216a05f60SBarry Smith    Limited by arguments such that `n` choose `k` can be represented by `PetscInt`
49316a05f60SBarry Smith 
49416a05f60SBarry Smith .seealso: `PetscDTEnumSubset()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()`
495fad4db65SToby Isaac M*/
PetscDTSubsetIndex(PetscInt n,PetscInt k,const PetscInt * subset,PetscInt * index)496d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
497d71ae5a4SJacob Faibussowitsch {
4985f80ce2aSJacob Faibussowitsch   PetscInt j = 0, Nk;
4991a989b97SToby Isaac 
50028222859SToby Isaac   PetscFunctionBegin;
50128222859SToby Isaac   *index = -1;
5029566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
5035f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
5041a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
5051a989b97SToby Isaac     PetscInt Nminusk      = Nk - Nminuskminus;
5061a989b97SToby Isaac 
5071a989b97SToby Isaac     if (subset[l] == i) {
5081a989b97SToby Isaac       l++;
5091a989b97SToby Isaac       Nk = Nminuskminus;
5101a989b97SToby Isaac     } else {
5111a989b97SToby Isaac       j += Nminuskminus;
5121a989b97SToby Isaac       Nk = Nminusk;
5131a989b97SToby Isaac     }
5141a989b97SToby Isaac   }
5151a989b97SToby Isaac   *index = j;
5163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5171a989b97SToby Isaac }
5181a989b97SToby Isaac 
519fad4db65SToby Isaac /*MC
52089521216SMatthew G. Knepley    PetscDTEnumSplit - 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.
521fad4db65SToby Isaac 
5224165533cSJose E. Roman    Input Parameters:
523fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
524fad4db65SToby Isaac .  k - an integer in [0, n]
525fad4db65SToby Isaac -  j - an index in [0, n choose k)
526fad4db65SToby Isaac 
5274165533cSJose E. Roman    Output Parameters:
528fad4db65SToby Isaac +  perm  - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set.
5292fe279fdSBarry Smith -  isOdd - if not `NULL`, return whether perm is an even or odd permutation.
530fad4db65SToby Isaac 
531fad4db65SToby Isaac    Level: beginner
532fad4db65SToby Isaac 
53316a05f60SBarry Smith    Note:
53416a05f60SBarry Smith    Limited by arguments such that `n` choose `k` can be represented by `PetscInt`
53516a05f60SBarry Smith 
53616a05f60SBarry Smith .seealso: `PetscDTEnumSubset()`, `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`,
53716a05f60SBarry Smith           `PetscDTPermIndex()`
538fad4db65SToby Isaac M*/
PetscDTEnumSplit(PetscInt n,PetscInt k,PetscInt j,PetscInt * perm,PeOp PetscBool * isOdd)539ce78bad3SBarry Smith static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PeOp PetscBool *isOdd)
540d71ae5a4SJacob Faibussowitsch {
5415f80ce2aSJacob Faibussowitsch   PetscInt  i, l, m, Nk, odd = 0;
5428e3a54c0SPierre Jolivet   PetscInt *subcomp = PetscSafePointerPlusOffset(perm, k);
5431a989b97SToby Isaac 
54428222859SToby Isaac   PetscFunctionBegin;
54528222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
5469566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
5471a989b97SToby Isaac   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
5481a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
5491a989b97SToby Isaac     PetscInt Nminusk      = Nk - Nminuskminus;
5501a989b97SToby Isaac 
5511a989b97SToby Isaac     if (j < Nminuskminus) {
552fad4db65SToby Isaac       perm[l++] = i;
5531a989b97SToby Isaac       Nk        = Nminuskminus;
5541a989b97SToby Isaac     } else {
5551a989b97SToby Isaac       subcomp[m++] = i;
5561a989b97SToby Isaac       j -= Nminuskminus;
5571a989b97SToby Isaac       odd ^= ((k - l) & 1);
5581a989b97SToby Isaac       Nk = Nminusk;
5591a989b97SToby Isaac     }
5601a989b97SToby Isaac   }
5615f80ce2aSJacob Faibussowitsch   for (; i < n; i++) subcomp[m++] = i;
5621a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
5633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5641a989b97SToby Isaac }
5651a989b97SToby Isaac 
566ce78bad3SBarry Smith struct _n_PetscTabulation {
567a5b23f4aSJose E. Roman   PetscInt    K;    /* Indicates a k-jet, namely tabulated derivatives up to order k */
56819815104SMartin Diehl   PetscInt    Nr;   /* The number of tabulation replicas (often 1) */
569ef0bb6c7SMatthew G. Knepley   PetscInt    Np;   /* The number of tabulation points in a replica */
570ef0bb6c7SMatthew G. Knepley   PetscInt    Nb;   /* The number of functions tabulated */
571ef0bb6c7SMatthew G. Knepley   PetscInt    Nc;   /* The number of function components */
572ef0bb6c7SMatthew G. Knepley   PetscInt    cdim; /* The coordinate dimension */
573ef0bb6c7SMatthew G. Knepley   PetscReal **T;    /* The tabulation T[K] of functions and their derivatives
574ef0bb6c7SMatthew G. Knepley                        T[0] = B[Nr*Np][Nb][Nc]:             The basis function values at quadrature points
575ef0bb6c7SMatthew G. Knepley                        T[1] = D[Nr*Np][Nb][Nc][cdim]:       The basis function derivatives at quadrature points
576ef0bb6c7SMatthew G. Knepley                        T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */
577ef0bb6c7SMatthew G. Knepley };
578ce78bad3SBarry Smith 
579ce78bad3SBarry Smith /*S
580ce78bad3SBarry Smith    PetscTabulation - PETSc object that manages tabulations for finite element methods.
581ce78bad3SBarry Smith 
582ce78bad3SBarry Smith    Level: intermediate
583ce78bad3SBarry Smith 
584ce78bad3SBarry Smith    Note:
585ce78bad3SBarry Smith    This is a pointer to a C struct, hence the data in it may be accessed directly.
586ce78bad3SBarry Smith 
587ce78bad3SBarry Smith    Fortran Note:
588ce78bad3SBarry Smith    Use `PetscTabulationGetData()` and `PetscTabulationRestoreData()` to access the arrays in the tabulation.
589ce78bad3SBarry Smith 
590ce78bad3SBarry Smith    Developer Note:
591ce78bad3SBarry Smith    TODO: put the meaning of the struct fields in this manual page
592ce78bad3SBarry Smith 
593ce78bad3SBarry Smith .seealso: `PetscTabulationDestroy()`, `PetscFECreateTabulation()`, `PetscFEGetCellTabulation()`
594ce78bad3SBarry Smith S*/
595ce78bad3SBarry Smith typedef struct _n_PetscTabulation *PetscTabulation;
596ef0bb6c7SMatthew G. Knepley 
597f8662bd6SBarry Smith /*S
598f8662bd6SBarry Smith   PetscProbFn - A prototype of a PDF or CDF used with PETSc probability operations whose names begin with `PetscProb` such as
599f8662bd6SBarry Smith   `PetscProbComputeKSStatistic()`.
600f8662bd6SBarry Smith 
601f8662bd6SBarry Smith   Calling Sequence:
602f8662bd6SBarry Smith + x      - input value
603f8662bd6SBarry Smith . scale  - scale factor, I don't know what this is for
604f8662bd6SBarry Smith - result - the value of the PDF or CDF at the input value
605f8662bd6SBarry Smith 
606f8662bd6SBarry Smith   Level: beginner
607f8662bd6SBarry Smith 
608f8662bd6SBarry Smith   Developer Note:
609f8662bd6SBarry Smith   Why does this take an array argument for `result` when it seems to be able to output a single value?
610f8662bd6SBarry Smith 
611f8662bd6SBarry Smith .seealso: `PetscProbComputeKSStatistic()`, `PetscProbComputeKSStatisticWeighted()`, `PetscPDFMaxwellBoltzmann1D()`
612f8662bd6SBarry Smith S*/
613f8662bd6SBarry Smith typedef PetscErrorCode PetscProbFn(const PetscReal x[], const PetscReal scale[], PetscReal result[]);
614f8662bd6SBarry Smith 
615f8662bd6SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscProbFn *PetscProbFunc PETSC_DEPRECATED_TYPEDEF(3, 24, 0, "PetscProbFn*", );
616d6685f55SMatthew G. Knepley 
6179371c9d4SSatish Balay typedef enum {
6189371c9d4SSatish Balay   DTPROB_DENSITY_CONSTANT,
6199371c9d4SSatish Balay   DTPROB_DENSITY_GAUSSIAN,
6209371c9d4SSatish Balay   DTPROB_DENSITY_MAXWELL_BOLTZMANN,
6219371c9d4SSatish Balay   DTPROB_NUM_DENSITY
6229371c9d4SSatish Balay } DTProbDensityType;
623d6685f55SMatthew G. Knepley PETSC_EXTERN const char *const DTProbDensityTypes[];
624d6685f55SMatthew G. Knepley 
625f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFMaxwellBoltzmann1D;
626f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscCDFMaxwellBoltzmann1D;
627f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFMaxwellBoltzmann2D;
628f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscCDFMaxwellBoltzmann2D;
629f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFMaxwellBoltzmann3D;
630f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscCDFMaxwellBoltzmann3D;
631f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFGaussian1D;
632f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscCDFGaussian1D;
633f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFSampleGaussian1D;
634f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFGaussian2D;
635f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFSampleGaussian2D;
636f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFGaussian3D;
637f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFSampleGaussian3D;
638f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFConstant1D;
639f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscCDFConstant1D;
640f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFSampleConstant1D;
641f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFConstant2D;
642f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscCDFConstant2D;
643f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFSampleConstant2D;
644f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFConstant3D;
645f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscCDFConstant3D;
646f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFSampleConstant3D;
647f8662bd6SBarry Smith PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFn **, PetscProbFn **, PetscProbFn **);
648d6685f55SMatthew G. Knepley 
649d6685f55SMatthew G. Knepley #include <petscvec.h>
650d6685f55SMatthew G. Knepley 
651f8662bd6SBarry Smith PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFn *, PetscReal *);
652f8662bd6SBarry Smith PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatisticWeighted(Vec, Vec, PetscProbFn *, PetscReal *);
653f8662bd6SBarry Smith PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatisticMagnitude(Vec, PetscProbFn *, PetscReal *);
654