DMPlexInsertBoundaryValuesEssentialField#
Insert boundary values into a local vector using a function of the coordinates and field data
Synopsis#
#include "petscdmplex.h"
PetscErrorCode DMPlexInsertBoundaryValuesEssentialField(DM dm, PetscReal time, Vec locU, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[], 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[], PetscInt, const PetscScalar[], PetscScalar[]), void *ctx, Vec locX)
Input Parameters#
dm - The
DM
, with aPetscDS
that matches the problem being constrainedtime - The time
locU - A local vector with the input solution values
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 componentslabel - The
DMLabel
defining constrained pointsnumids - The number of
DMLabel
ids for constrained pointsids - 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
, DMPlexInsertBoundaryValuesEssential()
, DMPlexInsertBoundaryValuesEssentialBdField()
, DMAddBoundary()
Level#
developer
Location#
Index of all DMPlex routines
Table of Contents for all manual pages
Index of all manual pages