DMPlexComputeL2DiffLocal#
This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
Synopsis#
#include "petscdmplex.h"
PetscErrorCode DMPlexComputeL2DiffLocal(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec localX, PetscReal *diff)
Collective
Input Parameters#
dm - The
DM
time - The time
funcs - The functions to evaluate for each field component
ctxs - Optional array of contexts to pass to each function, or
NULL
.localX - The coefficient vector u_h, a local vector
Output Parameter#
diff - The diff ||u - u_h||_2
See Also#
DMPlex: Unstructured Grids, DM
, DMPLEX
, DMProjectFunction()
, DMComputeL2FieldDiff()
, DMComputeL2GradientDiff()
Level#
developer
Location#
Index of all DMPlex routines
Table of Contents for all manual pages
Index of all manual pages