PETSc version 3.16.6
Fix/Edit manual page

PetscFVCreateTabulation

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

Synopsis

#include "petscfv.h" 
PetscErrorCode PetscFVCreateTabulation(PetscFV fvm, PetscInt nrepl, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation *T)
Not collective

Input Parameters

fvm - The PetscFV object
nrepl - The number of replicas
npoints - The number of tabulation points in a replica
points - The tabulation point coordinates
K - The order of derivative to tabulate

Output Parameter

T - The basis function values and derivative 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

PetscFECreateTabulation(), PetscTabulationDestroy(), PetscFEGetCellTabulation()

Level

intermediate

Location

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