petsc-3.9.4 2018-09-11
Report Typos and Errors

PetscFEGetTabulation

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

Synopsis

#include "petscfe.h" 
PetscErrorCode PetscFEGetTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H)
Not collective

Input Parameters

fem - The PetscFE object
npoints - The number of tabulation points
points - The tabulation point coordinates

Output Parameters

B - The basis function values at tabulation points
D - The basis function derivatives at tabulation points
H - The basis function second derivatives at tabulation 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

PetscFERestoreTabulation(), PetscFEGetDefaultTabulation()

Level

intermediate

Location

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