PetscFEGetFaceTabulation#

Returns the tabulation of the basis functions at the face quadrature points for each face of the reference cell

Synopsis#

#include "petscfe.h" 
PetscErrorCode PetscFEGetFaceTabulation(PetscFE fem, PetscInt k, PetscTabulation *Tf)

Not Collective

Input Parameters#

  • fem - The PetscFE object

  • k - The highest derivative we need to tabulate, very often 1

Output Parameter#

  • Tf - The basis function values and derivatives at face quadrature points

Note#

  T->T[0] = Bf[((f*Nq + q)*pdim + i)*Nc + c] is the value at point f,q for basis function i and component c
  T->T[1] = Df[(((f*Nq + q)*pdim + i)*Nc + c)*dim + d] is the derivative value at point f,q for basis function i, component c, in direction d
  T->T[2] = Hf[((((f*Nq + q)*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point f,q for basis function i, component c, in directions d and e

See Also#

PetscFE, PetscSpace, PetscDualSpace, PetscTabulation, PetscFEGetCellTabulation(), PetscFECreateTabulation(), 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