PetscDSUpdateBoundary#

Change a boundary condition for the model. The pointwise functions are used to provide boundary values for essential boundary conditions. In FEM, they are acting upon by dual basis functionals to generate FEM coefficients which are fixed. Natural boundary conditions signal to PETSc that boundary integrals should be performed, using the kernels from PetscDSSetBdResidual().

Synopsis#

#include "petscds.h" 
PetscErrorCode PetscDSUpdateBoundary(PetscDS ds, PetscInt bd, DMBoundaryConditionType type, const char name[], DMLabel label, PetscInt Nv, const PetscInt values[], PetscInt field, PetscInt Nc, const PetscInt comps[], void (*bcFunc)(void), void (*bcFunc_t)(void), void *ctx)

Input Parameters#

  • ds - The PetscDS object

  • bd - The boundary condition number

  • type - The type of condition, e.g. DM_BC_ESSENTIAL/DM_BC_ESSENTIAL_FIELD (Dirichlet), or DM_BC_NATURAL (Neumann)

  • name - The BC name

  • label - The label defining constrained points

  • Nv - The number of DMLabel ids for constrained points

  • values - An array of ids for constrained points

  • field - The field to constrain

  • Nc - The number of constrained field components

  • comps - An array of constrained component numbers

  • bcFunc - A pointwise function giving boundary values

  • bcFunc_t - A pointwise function giving the time derivative of the boundary values, or NULL

  • ctx - An optional user context for bcFunc

Note#

The boundary condition number is the order in which it was registered. The user can get the number of boundary conditions from PetscDSGetNumBoundary(). See PetscDSAddBoundary() for a description of the calling sequences for the callbacks.

See Also#

PetscDS, PetscWeakForm, DMBoundaryConditionType, PetscDSAddBoundary(), PetscDSGetBoundary(), PetscDSGetNumBoundary(), DMLabel

Level#

developer

Location#

src/dm/dt/interface/dtds.c


Edit on GitLab

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