petsc-3.14.6 2021-03-30
Report Typos and Errors

PetscDTJacobiEvalJet

Evaluate the jet (function and derivatives) of the Jacobi polynomials basis up to a given degree. The Jacobi polynomials with indices $\alpha$ and $\beta$ are orthogonal with respect to the weighted inner product $\langle f, g \rangle = \int_{-1}^1 (1+x)^{\alpha} (1-x)^{\beta) f(x) g(x) dx$.

Synopsis

#include "petscdt.h" 
PetscErrorCode PetscDTJacobiEvalJet(PetscReal alpha, PetscReal beta, PetscInt npoints, const PetscReal points[], PetscInt degree, PetscInt k, PetscReal p[])

Input Arguments

alpha - the left exponent of the weight
beta - the right exponetn of the weight
npoints - the number of points to evaluate the polynomials at
points - [npoints] array of point coordinates
degree - the maximm degree polynomial space to evaluate, (degree + 1) will be evaluated total.
k - the maximum derivative to evaluate in the jet, (k + 1) will be evaluated total.

Output Argments

p - an array containing the evaluations of the Jacobi polynomials's jets on the points. the size is (degree + 1) x (k + 1) x npoints, which also describes the order of the dimensions of this three-dimensional array: the first (slowest varying) dimension is polynomial degree; the second dimension is derivative order; the third (fastest varying) dimension is the index of the evaluation point.

See Also

PetscDTJacobiEval(), PetscDTPKDEvalJet()

Level

advanced

Location

src/dm/dt/interface/dt.c
Index of all DT routines
Table of Contents for all manual pages
Index of all manual pages