PetscDTPTrimmedEvalJet#

Evaluate the jet (function and derivatives) of a basis of the trimmed polynomial k-forms up to a given degree.

Synopsis#

#include "petscdt.h" 
PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt degree, PetscInt formDegree, PetscInt jetDegree, 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 trimmed polynomial space to evaluate. There are ((dim + degree) choose (dim + formDegree)) x ((degree + formDegree - 1) choose (formDegree)) polynomials in this space. (You can use PetscDTPTrimmedSize() to compute this size.)

  • formDegree - the degree of the form

  • jetDegree - the maximum order partial derivative to evaluate in the jet. There are ((dim + jetDegree) choose dim) partial derivatives in the jet. Choosing jetDegree = 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 PetscDTPTrimmedSize() x ((dim + formDegree) choose dim) x ((dim + k) choose dim) x npoints, which also describes the order of the dimensions of this four-dimensional array: the first (slowest varying) dimension is basis function index; the second dimension is component of the form; the third dimension is jet index; the fourth (fastest varying) dimension is the index of the evaluation point.

Notes#

The ordering of the basis functions is not graded, so the basis functions are not nested by degree like PetscDTPKDEvalJet(). The basis functions are not an L2-orthonormal basis on any particular domain.

The implementation is based on the description of the trimmed polynomials up to degree r as the direct sum of polynomials up to degree (r-1) and the Koszul differential applied to homogeneous polynomials of degree (r-1).

See Also#

PetscDTPKDEvalJet(), PetscDTPTrimmedSize()

Level#

advanced

Location#

src/dm/dt/interface/dt.c


Edit on GitLab

Index of all DT routines
Table of Contents for all manual pages
Index of all manual pages