PETSc version 3.16.6
Fix/Edit manual page

DMComputeL2Diff

This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.

Synopsis

#include "petscdm.h"          
#include "petscdmlabel.h"     
#include "petscds.h"     
PetscErrorCode DMComputeL2Diff(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)

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.
X - The coefficient vector u_h, a global vector

Output Parameter

diff - The diff ||u - u_h||_2

See Also

DMProjectFunction(), DMComputeL2FieldDiff(), DMComputeL2GradientDiff()

Level

developer

Location

src/dm/interface/dm.c

Examples

src/snes/tutorials/ex12.c.html
src/ts/tutorials/ex48.c.html
src/tao/tutorials/ex1.c.html
src/tao/tutorials/ex2.c.html

Implementations

DMComputeL2Diff_pforest(DM dm, PetscReal time, PetscErrorCode (**funcs) in src/dm/impls/forest/p4est/pforest.c
DMComputeL2Diff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs) in src/dm/impls/plex/plexfem.c

Index of all DM routines
Table of Contents for all manual pages
Index of all manual pages