PetscDTJacobiEval#

evaluate Jacobi polynomials for the weight function \((1.+x)^{\alpha} (1.-x)^{\beta}\) at a set of points at points

Synopsis#

#include "petscdt.h" 
PetscErrorCode PetscDTJacobiEval(PetscInt npoints, PetscReal alpha, PetscReal beta, const PetscReal *points, PetscInt ndegree, const PetscInt *degrees, PetscReal *B, PetscReal *D, PetscReal *D2)

Not Collective

Input Parameters#

  • npoints - number of spatial points to evaluate at

  • alpha - the left exponent > -1

  • beta - the right exponent > -1

  • points - array of locations to evaluate at

  • ndegree - number of basis degrees to evaluate

  • degrees - sorted array of degrees to evaluate

Output Parameters#

  • B - row-oriented basis evaluation matrix B[pointndegree + degree] (dimension npointsndegrees, allocated by caller) (or NULL)

  • D - row-oriented derivative evaluation matrix (or NULL)

  • D2 - row-oriented second derivative evaluation matrix (or NULL)

See Also#

PetscDTGaussQuadrature(), PetscDTLegendreEval()

Level#

intermediate

Location#

src/dm/dt/interface/dt.c


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