petsc-3.10.5 2019-03-28
DMPlexComputeL2DiffVec
This function computes the cellwise L_2 difference between a function u and an FEM interpolant solution u_h, and stores it in a Vec.
Synopsis
#include "petscdmplex.h"
PetscErrorCode DMPlexComputeL2DiffVec(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, Vec D)
Input Parameters
| dm | - The DM
|
| time | - The time
|
| funcs | - The functions to evaluate for each field component: NULL means that component does not contribute to error calculation
|
| ctxs | - Optional array of contexts to pass to each function, or NULL.
|
| X | - The coefficient vector u_h
|
Output Parameter
D -A Vec which holds the difference ||u - u_h||_2 for each cell
See Also
DMProjectFunction(), DMComputeL2Diff(), DMPlexComputeL2FieldDiff(), DMComputeL2GradientDiff()
Level
developer
Location
src/dm/impls/plex/plexfem.c
Index of all DMPLEX routines
Table of Contents for all manual pages
Index of all manual pages