:orphan: # PetscDTPKDEvalJet Evaluate the jet (function and derivatives) of the Proriol-Koornwinder-Dubiner (PKD) basis for the space of polynomials up to a given degree. The PKD basis is L2-orthonormal on the biunit simplex (which is used as the reference element for finite elements in PETSc), which makes it a stable basis to use for evaluating polynomials in that domain. ## Synopsis ``` #include "petscdt.h" PetscErrorCode PetscDTPKDEvalJet(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt degree, PetscInt k, PetscReal p[]) ``` ## Input Parameters - ***dim -*** the number of variables in the multivariate polynomials - ***npoints -*** the number of points to evaluate the polynomials at - ***points -*** [npoints x dim] array of point coordinates - ***degree -*** the degree (sum of degrees on the variables in a monomial) of the polynomial space to evaluate. There are ((dim + degree) choose dim) polynomials in this space. - ***k -*** the maximum order partial derivative to evaluate in the jet. There are (dim + k choose dim) partial derivatives in the jet. Choosing k = 0 means to evaluate just the function and no derivatives ## Output Parameter - ***p -*** an array containing the evaluations of the PKD polynomials' jets on the points. The size is ((dim + degree) choose dim) x ((dim + k) choose dim) x npoints, which also describes the order of the dimensions of this three-dimensional array: the first (slowest varying) dimension is basis function index; the second dimension is jet index; the third (fastest varying) dimension is the index of the evaluation point. ## Notes The ordering of the basis functions, and the ordering of the derivatives in the jet, both follow the graded ordering of `PetscDTIndexToGradedOrder()` and `PetscDTGradedOrderToIndex()`. For example, in 3D, the polynomial with leading monomial x^2,y^0,z^1, which has degree tuple (2,0,1), which by `PetscDTGradedOrderToIndex()` has index 12 (it is the 13th basis function in the space); the partial derivative $\partial_x \partial_z$ has order tuple (1,0,1), appears at index 6 in the jet (it is the 7th partial derivative in the jet). The implementation uses Kirby's singularity-free evaluation algorithm, https://doi.org/10.1145/1644001.1644006. ## See Also `PetscDTGradedOrderToIndex()`, `PetscDTIndexToGradedOrder()`, `PetscDTJacobiEvalJet()` ## Level advanced ## Location src/dm/dt/interface/dt.c --- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/dm/dt/interface/dt.c) [Index of all DT routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)