DMComputeL2GradientDiff#

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

Synopsis#

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

Collective

Input Parameters#

  • dm - The DM

  • time - The time

  • funcs - The gradient 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

  • n - The vector to project along

Output Parameter#

  • diff - The diff ||(grad u - grad u_h) . n||_2

Developer Notes#

This API is specific to only particular usage of DM

The notes need to provide some information about what has to be provided to the DM to be able to perform the computation.

See Also#

DM Basics, DM, DMProjectFunction(), DMComputeL2Diff(), DMComputeL2FieldDiff()

Level#

developer

Location#

src/dm/interface/dm.c

Implementations#

DMComputeL2GradientDiff_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