petsc-3.12.4 2020-02-04
PetscFVGetTabulation
Tabulates the basis functions, and perhaps derivatives, at the points provided.
Synopsis
#include "petscfv.h"
PetscErrorCode PetscFVGetTabulation(PetscFV fvm, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H)
Not collective
Input Parameters
| fvm | - The PetscFV 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
PetscFEGetTabulation(), PetscFERestoreTabulation(), PetscFEGetDefaultTabulation()
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