PetscFEComputeTabulation#

Tabulates the basis functions, and perhaps derivatives, at the points provided.

Synopsis#

#include "petscfe.h" 
PetscErrorCode PetscFEComputeTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T)

Not Collective

Input Parameters#

  • fem - The PetscFE object

  • npoints - The number of tabulation points

  • points - The tabulation point coordinates

  • K - The number of derivatives calculated

  • T - An existing tabulation object with enough allocated space

Output Parameter#

  • T - The basis function values and derivatives at tabulation points

Note#

  T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
  T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
  T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis function i, component c, in directions d and e

See Also#

PetscTabulation, PetscFEGetCellTabulation(), PetscTabulationDestroy()

Level#

intermediate

Location#

src/dm/dt/fe/interface/fe.c


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