petsc-3.10.5 2019-03-28
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)
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
See Also
DMProjectFunction(), DMComputeL2Diff()
Level
developer
Location
src/dm/interface/dm.c
Implementations
DMComputeL2GradientDiff_DA(DM dm, PetscReal time, PetscErrorCode (**funcs) in src/dm/impls/da/dalocal.c
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