petsc-3.11.4 2019-09-28
Report Typos and Errors

PetscFEGetDefaultTabulation

Returns the tabulation of the basis functions at the quadrature points

Synopsis

#include "petscfe.h" 
PetscErrorCode PetscFEGetDefaultTabulation(PetscFE fem, PetscReal **B, PetscReal **D, PetscReal **H)
Not collective

Input Parameter

fem -The PetscFE object

Output Parameters

B - The basis function values at quadrature points
D - The basis function derivatives at quadrature points
H - The basis function second derivatives at quadrature points

Note

B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
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

PetscFEGetTabulation(), PetscFERestoreTabulation()

Level

intermediate

Location

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