xref: /petsc/include/petscdt.h (revision f8662bd6eba1da879f288a6933281c3de08ca1fb)
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 *);
13637045ce4SJed Brown 
137b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
138d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
139d6685f55SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
140b3c0f97bSTom Klotz 
141916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *);
142916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
143916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
144916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
145916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
146916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
147916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
148916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
149916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
150916e780bShannah_mairs 
1512dce792eSToby Isaac /*MC
1522dce792eSToby Isaac   PETSC_FORM_DEGREE_UNDEFINED - Indicates that a field does not have
1532dce792eSToby Isaac   a well-defined form degree in exterior calculus.
1542dce792eSToby Isaac 
1552dce792eSToby Isaac   Level: advanced
1562dce792eSToby Isaac 
1572dce792eSToby Isaac .seealso: `PetscDTAltV`, `PetscDualSpaceGetFormDegree()`
1582dce792eSToby Isaac M*/
1592dce792eSToby Isaac #define PETSC_FORM_DEGREE_UNDEFINED PETSC_INT_MIN
1602dce792eSToby Isaac 
1611a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1621a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1631a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
1641a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
1651a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
1661a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
1671a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
168dda711d0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
1691a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
1701a989b97SToby Isaac 
171d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt, PetscInt, const PetscInt[], PetscInt *);
172d4afb720SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt, PetscInt, PetscInt, PetscInt[]);
173fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTGradedOrderToIndex(PetscInt, const PetscInt[], PetscInt *);
174fbdc3dfeSToby Isaac PETSC_EXTERN PetscErrorCode PetscDTIndexToGradedOrder(PetscInt, PetscInt, PetscInt[]);
175d4afb720SToby Isaac 
176fad4db65SToby Isaac #if defined(PETSC_USE_64BIT_INDICES)
177fad4db65SToby Isaac   #define PETSC_FACTORIAL_MAX 20
178fad4db65SToby Isaac   #define PETSC_BINOMIAL_MAX  61
179fad4db65SToby Isaac #else
180fad4db65SToby Isaac   #define PETSC_FACTORIAL_MAX 12
181fad4db65SToby Isaac   #define PETSC_BINOMIAL_MAX  29
182fad4db65SToby Isaac #endif
183fad4db65SToby Isaac 
184fad4db65SToby Isaac /*MC
185fad4db65SToby Isaac    PetscDTFactorial - Approximate n! as a real number
186fad4db65SToby Isaac 
1874165533cSJose E. Roman    Input Parameter:
188fad4db65SToby Isaac .  n - a non-negative integer
189fad4db65SToby Isaac 
1904165533cSJose E. Roman    Output Parameter:
191fad4db65SToby Isaac .  factorial - n!
192fad4db65SToby Isaac 
193fad4db65SToby Isaac    Level: beginner
19416a05f60SBarry Smith 
19516a05f60SBarry Smith .seealso: `PetscDTFactorialInt()`, `PetscDTBinomialInt()`, `PetscDTBinomial()`
196fad4db65SToby Isaac M*/
197d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial)
198d71ae5a4SJacob Faibussowitsch {
199fad4db65SToby Isaac   PetscReal f = 1.0;
200fad4db65SToby Isaac 
201fad4db65SToby Isaac   PetscFunctionBegin;
202e2ab39ccSLisandro Dalcin   *factorial = -1.0;
20363a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %" PetscInt_FMT, n);
2045f80ce2aSJacob Faibussowitsch   for (PetscInt i = 1; i < n + 1; ++i) f *= (PetscReal)i;
205fad4db65SToby Isaac   *factorial = f;
2063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
207fad4db65SToby Isaac }
208fad4db65SToby Isaac 
209fad4db65SToby Isaac /*MC
210fad4db65SToby Isaac    PetscDTFactorialInt - Compute n! as an integer
211fad4db65SToby Isaac 
2124165533cSJose E. Roman    Input Parameter:
213fad4db65SToby Isaac .  n - a non-negative integer
214fad4db65SToby Isaac 
2154165533cSJose E. Roman    Output Parameter:
216fad4db65SToby Isaac .  factorial - n!
217fad4db65SToby Isaac 
218fad4db65SToby Isaac    Level: beginner
219fad4db65SToby Isaac 
22016a05f60SBarry Smith    Note:
22195bd0b28SBarry 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.
22216a05f60SBarry Smith 
22316a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTBinomialInt()`, `PetscDTBinomial()`
224fad4db65SToby Isaac M*/
225d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial)
226d71ae5a4SJacob Faibussowitsch {
227fad4db65SToby Isaac   PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600};
228fad4db65SToby Isaac 
22928222859SToby Isaac   PetscFunctionBegin;
23028222859SToby Isaac   *factorial = -1;
23163a3b9bcSJacob 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);
232fad4db65SToby Isaac   if (n <= 12) {
233fad4db65SToby Isaac     *factorial = facLookup[n];
234fad4db65SToby Isaac   } else {
235fad4db65SToby Isaac     PetscInt f = facLookup[12];
236fad4db65SToby Isaac     PetscInt i;
237fad4db65SToby Isaac 
238fad4db65SToby Isaac     for (i = 13; i < n + 1; ++i) f *= i;
239fad4db65SToby Isaac     *factorial = f;
240fad4db65SToby Isaac   }
2413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
242fad4db65SToby Isaac }
243fad4db65SToby Isaac 
244fad4db65SToby Isaac /*MC
24595bd0b28SBarry Smith    PetscDTBinomial - Approximate the binomial coefficient `n` choose `k`
246fad4db65SToby Isaac 
2474165533cSJose E. Roman    Input Parameters:
248fad4db65SToby Isaac +  n - a non-negative integer
24995bd0b28SBarry Smith -  k - an integer between 0 and `n`, inclusive
250fad4db65SToby Isaac 
2514165533cSJose E. Roman    Output Parameter:
25295bd0b28SBarry Smith .  binomial - approximation of the binomial coefficient `n` choose `k`
253fad4db65SToby Isaac 
254fad4db65SToby Isaac    Level: beginner
25516a05f60SBarry Smith 
25616a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`
257fad4db65SToby Isaac M*/
258d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
259d71ae5a4SJacob Faibussowitsch {
2601a989b97SToby Isaac   PetscFunctionBeginHot;
261e2ab39ccSLisandro Dalcin   *binomial = -1.0;
26263a3b9bcSJacob 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);
2631a989b97SToby Isaac   if (n <= 3) {
2649371c9d4SSatish Balay     PetscInt binomLookup[4][4] = {
2659371c9d4SSatish Balay       {1, 0, 0, 0},
2669371c9d4SSatish Balay       {1, 1, 0, 0},
2679371c9d4SSatish Balay       {1, 2, 1, 0},
2689371c9d4SSatish Balay       {1, 3, 3, 1}
2699371c9d4SSatish Balay     };
2701a989b97SToby Isaac 
271e2ab39ccSLisandro Dalcin     *binomial = (PetscReal)binomLookup[n][k];
2721a989b97SToby Isaac   } else {
273e2ab39ccSLisandro Dalcin     PetscReal binom = 1.0;
2741a989b97SToby Isaac 
2751a989b97SToby Isaac     k = PetscMin(k, n - k);
2765f80ce2aSJacob Faibussowitsch     for (PetscInt i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1);
2771a989b97SToby Isaac     *binomial = binom;
2781a989b97SToby Isaac   }
2793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2801a989b97SToby Isaac }
2811a989b97SToby Isaac 
282fad4db65SToby Isaac /*MC
28395bd0b28SBarry Smith    PetscDTBinomialInt - Compute the binomial coefficient `n` choose `k`
284fad4db65SToby Isaac 
28597bb3fdcSJose E. Roman    Input Parameters:
286fad4db65SToby Isaac +  n - a non-negative integer
28795bd0b28SBarry Smith -  k - an integer between 0 and `n`, inclusive
288fad4db65SToby Isaac 
28997bb3fdcSJose E. Roman    Output Parameter:
29095bd0b28SBarry Smith .  binomial - the binomial coefficient `n` choose `k`
291fad4db65SToby Isaac 
292fad4db65SToby Isaac    Level: beginner
29316a05f60SBarry Smith 
29416a05f60SBarry Smith    Note:
29516a05f60SBarry Smith    This is limited by integers that can be represented by `PetscInt`.
29616a05f60SBarry Smith 
29716a05f60SBarry Smith    Use `PetscDTBinomial()` for real number approximations of larger values
29816a05f60SBarry Smith 
29916a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTEnumPerm()`
300fad4db65SToby Isaac M*/
301d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial)
302d71ae5a4SJacob Faibussowitsch {
30328222859SToby Isaac   PetscInt bin;
30428222859SToby Isaac 
30528222859SToby Isaac   PetscFunctionBegin;
30628222859SToby Isaac   *binomial = -1;
30763a3b9bcSJacob 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);
30863a3b9bcSJacob 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);
309fad4db65SToby Isaac   if (n <= 3) {
3109371c9d4SSatish Balay     PetscInt binomLookup[4][4] = {
3119371c9d4SSatish Balay       {1, 0, 0, 0},
3129371c9d4SSatish Balay       {1, 1, 0, 0},
3139371c9d4SSatish Balay       {1, 2, 1, 0},
3149371c9d4SSatish Balay       {1, 3, 3, 1}
3159371c9d4SSatish Balay     };
316fad4db65SToby Isaac 
31728222859SToby Isaac     bin = binomLookup[n][k];
318fad4db65SToby Isaac   } else {
319fad4db65SToby Isaac     PetscInt binom = 1;
320fad4db65SToby Isaac 
321fad4db65SToby Isaac     k = PetscMin(k, n - k);
3225f80ce2aSJacob Faibussowitsch     for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
32328222859SToby Isaac     bin = binom;
324fad4db65SToby Isaac   }
32528222859SToby Isaac   *binomial = bin;
3263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
327fad4db65SToby Isaac }
328fad4db65SToby Isaac 
329ce78bad3SBarry Smith /* the following inline routines should be not be inline routines and then Fortran binding can be built automatically */
330ce78bad3SBarry Smith #define PeOp
331ce78bad3SBarry Smith 
332fad4db65SToby Isaac /*MC
33316a05f60SBarry Smith    PetscDTEnumPerm - Get a permutation of `n` integers from its encoding into the integers [0, n!) as a sequence of swaps.
334fad4db65SToby Isaac 
3354165533cSJose E. Roman    Input Parameters:
336fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
3378cd1e013SToby Isaac -  k - an integer in [0, n!)
338fad4db65SToby Isaac 
3394165533cSJose E. Roman    Output Parameters:
340fad4db65SToby Isaac +  perm  - the permuted list of the integers [0, ..., n-1]
3412fe279fdSBarry Smith -  isOdd - if not `NULL`, returns whether the permutation used an even or odd number of swaps.
342fad4db65SToby Isaac 
34316a05f60SBarry Smith    Level: intermediate
344fad4db65SToby Isaac 
34516a05f60SBarry Smith    Notes:
34616a05f60SBarry Smith    A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation,
34716a05f60SBarry Smith    by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in
34816a05f60SBarry Smith    some position j >= i.  This swap is encoded as the difference (j - i).  The difference d_i at step i is less than
34916a05f60SBarry Smith    (n - i).  This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number
35016a05f60SBarry Smith    (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}.
35116a05f60SBarry Smith 
35216a05f60SBarry 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.
35316a05f60SBarry Smith 
35416a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTPermIndex()`
355fad4db65SToby Isaac M*/
356ce78bad3SBarry Smith static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PeOp PetscBool *isOdd)
357d71ae5a4SJacob Faibussowitsch {
3581a989b97SToby Isaac   PetscInt  odd = 0;
3591a989b97SToby Isaac   PetscInt  i;
360fad4db65SToby Isaac   PetscInt  work[PETSC_FACTORIAL_MAX];
361fad4db65SToby Isaac   PetscInt *w;
3621a989b97SToby Isaac 
36328222859SToby Isaac   PetscFunctionBegin;
36428222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
36563a3b9bcSJacob 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);
366810441c8SPierre Jolivet   if (n >= 2) {
367fad4db65SToby Isaac     w = &work[n - 2];
3681a989b97SToby Isaac     for (i = 2; i <= n; i++) {
3691a989b97SToby Isaac       *(w--) = k % i;
3701a989b97SToby Isaac       k /= i;
3711a989b97SToby Isaac     }
372810441c8SPierre Jolivet   }
3731a989b97SToby Isaac   for (i = 0; i < n; i++) perm[i] = i;
3741a989b97SToby Isaac   for (i = 0; i < n - 1; i++) {
3751a989b97SToby Isaac     PetscInt s    = work[i];
3761a989b97SToby Isaac     PetscInt swap = perm[i];
3771a989b97SToby Isaac 
3781a989b97SToby Isaac     perm[i]     = perm[i + s];
3791a989b97SToby Isaac     perm[i + s] = swap;
3801a989b97SToby Isaac     odd ^= (!!s);
3811a989b97SToby Isaac   }
3821a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
3833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3841a989b97SToby Isaac }
3851a989b97SToby Isaac 
386fad4db65SToby Isaac /*MC
38716a05f60SBarry Smith    PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!).  This inverts `PetscDTEnumPerm()`.
3888cd1e013SToby Isaac 
3894165533cSJose E. Roman    Input Parameters:
3908cd1e013SToby Isaac +  n    - a non-negative integer (see note about limits below)
3918cd1e013SToby Isaac -  perm - the permuted list of the integers [0, ..., n-1]
3928cd1e013SToby Isaac 
3934165533cSJose E. Roman    Output Parameters:
3948cd1e013SToby Isaac +  k     - an integer in [0, n!)
3952fe279fdSBarry Smith -  isOdd - if not `NULL`, returns whether the permutation used an even or odd number of swaps.
3968cd1e013SToby Isaac 
3978cd1e013SToby Isaac    Level: beginner
39816a05f60SBarry Smith 
39916a05f60SBarry Smith    Note:
40016a05f60SBarry 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.
40116a05f60SBarry Smith 
40216a05f60SBarry Smith .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`
4038cd1e013SToby Isaac M*/
404ce78bad3SBarry Smith static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PeOp PetscBool *isOdd)
405d71ae5a4SJacob Faibussowitsch {
4068cd1e013SToby Isaac   PetscInt odd = 0;
4078cd1e013SToby Isaac   PetscInt i, idx;
4088cd1e013SToby Isaac   PetscInt work[PETSC_FACTORIAL_MAX];
4098cd1e013SToby Isaac   PetscInt iwork[PETSC_FACTORIAL_MAX];
4108cd1e013SToby Isaac 
4118cd1e013SToby Isaac   PetscFunctionBeginHot;
41228222859SToby Isaac   *k = -1;
41328222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
41463a3b9bcSJacob 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);
4158cd1e013SToby Isaac   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
4168cd1e013SToby Isaac   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
4178cd1e013SToby Isaac   for (idx = 0, i = 0; i < n - 1; i++) {
4188cd1e013SToby Isaac     PetscInt j    = perm[i];
4198cd1e013SToby Isaac     PetscInt icur = work[i];
4208cd1e013SToby Isaac     PetscInt jloc = iwork[j];
4218cd1e013SToby Isaac     PetscInt diff = jloc - i;
4228cd1e013SToby Isaac 
4238cd1e013SToby Isaac     idx = idx * (n - i) + diff;
4248cd1e013SToby Isaac     /* swap (i, jloc) */
4258cd1e013SToby Isaac     work[i]     = j;
4268cd1e013SToby Isaac     work[jloc]  = icur;
4278cd1e013SToby Isaac     iwork[j]    = i;
4288cd1e013SToby Isaac     iwork[icur] = jloc;
4298cd1e013SToby Isaac     odd ^= (!!diff);
4308cd1e013SToby Isaac   }
4318cd1e013SToby Isaac   *k = idx;
4328cd1e013SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
4333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4348cd1e013SToby Isaac }
4358cd1e013SToby Isaac 
4368cd1e013SToby Isaac /*MC
437fad4db65SToby Isaac    PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k).
438fad4db65SToby Isaac    The encoding is in lexicographic order.
439fad4db65SToby Isaac 
4404165533cSJose E. Roman    Input Parameters:
441fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
442fad4db65SToby Isaac .  k - an integer in [0, n]
443fad4db65SToby Isaac -  j - an index in [0, n choose k)
444fad4db65SToby Isaac 
4454165533cSJose E. Roman    Output Parameter:
446fad4db65SToby Isaac .  subset - the jth subset of size k of the integers [0, ..., n - 1]
447fad4db65SToby Isaac 
448fad4db65SToby Isaac    Level: beginner
449fad4db65SToby Isaac 
45016a05f60SBarry Smith    Note:
45116a05f60SBarry Smith    Limited by arguments such that `n` choose `k` can be represented by `PetscInt`
45216a05f60SBarry Smith 
45316a05f60SBarry Smith .seealso: `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()`
454fad4db65SToby Isaac M*/
455d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
456d71ae5a4SJacob Faibussowitsch {
4575f80ce2aSJacob Faibussowitsch   PetscInt Nk;
4581a989b97SToby Isaac 
4591a989b97SToby Isaac   PetscFunctionBeginHot;
4609566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
4615f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
4621a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
4631a989b97SToby Isaac     PetscInt Nminusk      = Nk - Nminuskminus;
4641a989b97SToby Isaac 
4651a989b97SToby Isaac     if (j < Nminuskminus) {
4661a989b97SToby Isaac       subset[l++] = i;
4671a989b97SToby Isaac       Nk          = Nminuskminus;
4681a989b97SToby Isaac     } else {
4691a989b97SToby Isaac       j -= Nminuskminus;
4701a989b97SToby Isaac       Nk = Nminusk;
4711a989b97SToby Isaac     }
4721a989b97SToby Isaac   }
4733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4741a989b97SToby Isaac }
4751a989b97SToby Isaac 
476fad4db65SToby Isaac /*MC
47787497f52SBarry 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.
47887497f52SBarry Smith    This is the inverse of `PetscDTEnumSubset`.
479fad4db65SToby Isaac 
4804165533cSJose E. Roman    Input Parameters:
481fad4db65SToby Isaac +  n      - a non-negative integer (see note about limits below)
482fad4db65SToby Isaac .  k      - an integer in [0, n]
483fad4db65SToby Isaac -  subset - an ordered subset of the integers [0, ..., n - 1]
484fad4db65SToby Isaac 
4854165533cSJose E. Roman    Output Parameter:
486fad4db65SToby Isaac .  index - the rank of the subset in lexicographic order
487fad4db65SToby Isaac 
488fad4db65SToby Isaac    Level: beginner
489fad4db65SToby Isaac 
49016a05f60SBarry Smith    Note:
49116a05f60SBarry Smith    Limited by arguments such that `n` choose `k` can be represented by `PetscInt`
49216a05f60SBarry Smith 
49316a05f60SBarry Smith .seealso: `PetscDTEnumSubset()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()`
494fad4db65SToby Isaac M*/
495d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
496d71ae5a4SJacob Faibussowitsch {
4975f80ce2aSJacob Faibussowitsch   PetscInt j = 0, Nk;
4981a989b97SToby Isaac 
49928222859SToby Isaac   PetscFunctionBegin;
50028222859SToby Isaac   *index = -1;
5019566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
5025f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
5031a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
5041a989b97SToby Isaac     PetscInt Nminusk      = Nk - Nminuskminus;
5051a989b97SToby Isaac 
5061a989b97SToby Isaac     if (subset[l] == i) {
5071a989b97SToby Isaac       l++;
5081a989b97SToby Isaac       Nk = Nminuskminus;
5091a989b97SToby Isaac     } else {
5101a989b97SToby Isaac       j += Nminuskminus;
5111a989b97SToby Isaac       Nk = Nminusk;
5121a989b97SToby Isaac     }
5131a989b97SToby Isaac   }
5141a989b97SToby Isaac   *index = j;
5153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5161a989b97SToby Isaac }
5171a989b97SToby Isaac 
518fad4db65SToby Isaac /*MC
51989521216SMatthew 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.
520fad4db65SToby Isaac 
5214165533cSJose E. Roman    Input Parameters:
522fad4db65SToby Isaac +  n - a non-negative integer (see note about limits below)
523fad4db65SToby Isaac .  k - an integer in [0, n]
524fad4db65SToby Isaac -  j - an index in [0, n choose k)
525fad4db65SToby Isaac 
5264165533cSJose E. Roman    Output Parameters:
527fad4db65SToby Isaac +  perm  - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set.
5282fe279fdSBarry Smith -  isOdd - if not `NULL`, return whether perm is an even or odd permutation.
529fad4db65SToby Isaac 
530fad4db65SToby Isaac    Level: beginner
531fad4db65SToby Isaac 
53216a05f60SBarry Smith    Note:
53316a05f60SBarry Smith    Limited by arguments such that `n` choose `k` can be represented by `PetscInt`
53416a05f60SBarry Smith 
53516a05f60SBarry Smith .seealso: `PetscDTEnumSubset()`, `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`,
53616a05f60SBarry Smith           `PetscDTPermIndex()`
537fad4db65SToby Isaac M*/
538ce78bad3SBarry Smith static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PeOp PetscBool *isOdd)
539d71ae5a4SJacob Faibussowitsch {
5405f80ce2aSJacob Faibussowitsch   PetscInt  i, l, m, Nk, odd = 0;
5418e3a54c0SPierre Jolivet   PetscInt *subcomp = PetscSafePointerPlusOffset(perm, k);
5421a989b97SToby Isaac 
54328222859SToby Isaac   PetscFunctionBegin;
54428222859SToby Isaac   if (isOdd) *isOdd = PETSC_FALSE;
5459566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n, k, &Nk));
5461a989b97SToby Isaac   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
5471a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
5481a989b97SToby Isaac     PetscInt Nminusk      = Nk - Nminuskminus;
5491a989b97SToby Isaac 
5501a989b97SToby Isaac     if (j < Nminuskminus) {
551fad4db65SToby Isaac       perm[l++] = i;
5521a989b97SToby Isaac       Nk        = Nminuskminus;
5531a989b97SToby Isaac     } else {
5541a989b97SToby Isaac       subcomp[m++] = i;
5551a989b97SToby Isaac       j -= Nminuskminus;
5561a989b97SToby Isaac       odd ^= ((k - l) & 1);
5571a989b97SToby Isaac       Nk = Nminusk;
5581a989b97SToby Isaac     }
5591a989b97SToby Isaac   }
5605f80ce2aSJacob Faibussowitsch   for (; i < n; i++) subcomp[m++] = i;
5611a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
5623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5631a989b97SToby Isaac }
5641a989b97SToby Isaac 
565ce78bad3SBarry Smith struct _n_PetscTabulation {
566a5b23f4aSJose E. Roman   PetscInt    K;    /* Indicates a k-jet, namely tabulated derivatives up to order k */
56719815104SMartin Diehl   PetscInt    Nr;   /* The number of tabulation replicas (often 1) */
568ef0bb6c7SMatthew G. Knepley   PetscInt    Np;   /* The number of tabulation points in a replica */
569ef0bb6c7SMatthew G. Knepley   PetscInt    Nb;   /* The number of functions tabulated */
570ef0bb6c7SMatthew G. Knepley   PetscInt    Nc;   /* The number of function components */
571ef0bb6c7SMatthew G. Knepley   PetscInt    cdim; /* The coordinate dimension */
572ef0bb6c7SMatthew G. Knepley   PetscReal **T;    /* The tabulation T[K] of functions and their derivatives
573ef0bb6c7SMatthew G. Knepley                        T[0] = B[Nr*Np][Nb][Nc]:             The basis function values at quadrature points
574ef0bb6c7SMatthew G. Knepley                        T[1] = D[Nr*Np][Nb][Nc][cdim]:       The basis function derivatives at quadrature points
575ef0bb6c7SMatthew G. Knepley                        T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */
576ef0bb6c7SMatthew G. Knepley };
577ce78bad3SBarry Smith 
578ce78bad3SBarry Smith /*S
579ce78bad3SBarry Smith    PetscTabulation - PETSc object that manages tabulations for finite element methods.
580ce78bad3SBarry Smith 
581ce78bad3SBarry Smith    Level: intermediate
582ce78bad3SBarry Smith 
583ce78bad3SBarry Smith    Note:
584ce78bad3SBarry Smith    This is a pointer to a C struct, hence the data in it may be accessed directly.
585ce78bad3SBarry Smith 
586ce78bad3SBarry Smith    Fortran Note:
587ce78bad3SBarry Smith    Use `PetscTabulationGetData()` and `PetscTabulationRestoreData()` to access the arrays in the tabulation.
588ce78bad3SBarry Smith 
589ce78bad3SBarry Smith    Developer Note:
590ce78bad3SBarry Smith    TODO: put the meaning of the struct fields in this manual page
591ce78bad3SBarry Smith 
592ce78bad3SBarry Smith .seealso: `PetscTabulationDestroy()`, `PetscFECreateTabulation()`, `PetscFEGetCellTabulation()`
593ce78bad3SBarry Smith S*/
594ce78bad3SBarry Smith typedef struct _n_PetscTabulation *PetscTabulation;
595ef0bb6c7SMatthew G. Knepley 
596*f8662bd6SBarry Smith /*S
597*f8662bd6SBarry Smith   PetscProbFn - A prototype of a PDF or CDF used with PETSc probability operations whose names begin with `PetscProb` such as
598*f8662bd6SBarry Smith   `PetscProbComputeKSStatistic()`.
599*f8662bd6SBarry Smith 
600*f8662bd6SBarry Smith   Calling Sequence:
601*f8662bd6SBarry Smith + x      - input value
602*f8662bd6SBarry Smith . scale  - scale factor, I don't know what this is for
603*f8662bd6SBarry Smith - result - the value of the PDF or CDF at the input value
604*f8662bd6SBarry Smith 
605*f8662bd6SBarry Smith   Level: beginner
606*f8662bd6SBarry Smith 
607*f8662bd6SBarry Smith   Developer Note:
608*f8662bd6SBarry Smith   Why does this take an array argument for `result` when it seems to be able to output a single value?
609*f8662bd6SBarry Smith 
610*f8662bd6SBarry Smith .seealso: `PetscProbComputeKSStatistic()`, `PetscProbComputeKSStatisticWeighted()`, `PetscPDFMaxwellBoltzmann1D()`
611*f8662bd6SBarry Smith S*/
612*f8662bd6SBarry Smith typedef PetscErrorCode PetscProbFn(const PetscReal x[], const PetscReal scale[], PetscReal result[]);
613*f8662bd6SBarry Smith 
614*f8662bd6SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscProbFn *PetscProbFunc PETSC_DEPRECATED_TYPEDEF(3, 24, 0, "PetscProbFn*", );
615d6685f55SMatthew G. Knepley 
6169371c9d4SSatish Balay typedef enum {
6179371c9d4SSatish Balay   DTPROB_DENSITY_CONSTANT,
6189371c9d4SSatish Balay   DTPROB_DENSITY_GAUSSIAN,
6199371c9d4SSatish Balay   DTPROB_DENSITY_MAXWELL_BOLTZMANN,
6209371c9d4SSatish Balay   DTPROB_NUM_DENSITY
6219371c9d4SSatish Balay } DTProbDensityType;
622d6685f55SMatthew G. Knepley PETSC_EXTERN const char *const DTProbDensityTypes[];
623d6685f55SMatthew G. Knepley 
624*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFMaxwellBoltzmann1D;
625*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscCDFMaxwellBoltzmann1D;
626*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFMaxwellBoltzmann2D;
627*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscCDFMaxwellBoltzmann2D;
628*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFMaxwellBoltzmann3D;
629*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscCDFMaxwellBoltzmann3D;
630*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFGaussian1D;
631*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscCDFGaussian1D;
632*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFSampleGaussian1D;
633*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFGaussian2D;
634*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFSampleGaussian2D;
635*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFGaussian3D;
636*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFSampleGaussian3D;
637*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFConstant1D;
638*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscCDFConstant1D;
639*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFSampleConstant1D;
640*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFConstant2D;
641*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscCDFConstant2D;
642*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFSampleConstant2D;
643*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFConstant3D;
644*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscCDFConstant3D;
645*f8662bd6SBarry Smith PETSC_EXTERN PetscProbFn    PetscPDFSampleConstant3D;
646*f8662bd6SBarry Smith PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFn **, PetscProbFn **, PetscProbFn **);
647d6685f55SMatthew G. Knepley 
648d6685f55SMatthew G. Knepley #include <petscvec.h>
649d6685f55SMatthew G. Knepley 
650*f8662bd6SBarry Smith PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFn *, PetscReal *);
651*f8662bd6SBarry Smith PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatisticWeighted(Vec, Vec, PetscProbFn *, PetscReal *);
652*f8662bd6SBarry Smith PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatisticMagnitude(Vec, PetscProbFn *, PetscReal *);
653