DMPlexComputeBdIntegral#

Form the integral over the specified boundary from the global input X using pointwise functions specified by the user

Synopsis#

#include "petscdmplex.h"   
PetscErrorCode DMPlexComputeBdIntegral(DM dm, Vec X, DMLabel label, PetscInt numVals, const PetscInt vals[], void (*func)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscScalar *integral, void *user)

Input Parameters#

  • dm - The mesh

  • X - Global input vector

  • label - The boundary DMLabel

  • numVals - The number of label values to use, or PETSC_DETERMINE for all values

  • vals - The label values to use, or NULL for all values

  • func - The function to integrate along the boundary

  • user - The user context

Output Parameter#

  • integral - Integral for each field

See Also#

DMPlex: Unstructured Grids, DM, DMPLEX, DMPlexComputeIntegralFEM(), DMPlexComputeBdResidualFEM()

Level#

developer

Location#

src/dm/impls/plex/plexfem.c

Examples#

src/snes/tutorials/ex12.c

Implementations#

DMPlexComputeBdIntegral_Internal(DM dm, Vec locX, IS pointIS, void (*func) in src/dm/impls/plex/plexfem.c


Edit on GitLab

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