DMPlexInsertBoundaryValuesEssential#

Insert boundary values into a local vector using a function of the coordinates

Synopsis#

#include "petscdmplex.h"   
PetscErrorCode DMPlexInsertBoundaryValuesEssential(DM dm, PetscReal time, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *ctx, Vec locX)

Input Parameters#

  • dm - The DM, with a PetscDS that matches the problem being constrained

  • time - The time

  • field - The field to constrain

  • Nc - The number of constrained field components, or 0 for all components

  • comps - An array of constrained component numbers, or NULL for all components

  • label - The DMLabel defining constrained points

  • numids - The number of DMLabel ids for constrained points

  • ids - An array of ids for constrained points

  • func - A pointwise function giving boundary values

  • ctx - An optional user context for bcFunc

Output Parameter#

  • locX - A local vector to receives the boundary values

See Also#

DMPlex: Unstructured Grids, DM, DMPLEX, DMLabel, DMPlexInsertBoundaryValuesEssentialField(), DMPlexInsertBoundaryValuesEssentialBdField(), DMAddBoundary()

Level#

developer

Location#

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