PetscDSAddBoundaryByName#
Add a boundary condition to the model.
Synopsis#
#include "petscds.h" 
PetscErrorCode PetscDSAddBoundaryByName(PetscDS ds, DMBoundaryConditionType type, const char name[], const char lname[], PetscInt Nv, const PetscInt values[], PetscInt field, PetscInt Nc, const PetscInt comps[], void (*bcFunc)(void), void (*bcFunc_t)(void), void *ctx, PetscInt *bd)
Collective
Input Parameters#
- ds - The - PetscDSobject
- type - The type of condition, e.g. - DM_BC_ESSENTIAL/- DM_BC_ESSENTIAL_FIELD(Dirichlet), or- DM_BC_NATURAL(Neumann)
- name - The BC name 
- lname - The naem of the label defining constrained points 
- Nv - The number of - DMLabelvalues for constrained points
- values - An array of label values for constrained points 
- field - The field to constrain 
- Nc - The number of constrained field components (0 will constrain all fields) 
- 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 
Output Parameter#
- bd - The boundary number 
Options Database Keys#
- -bc_ - - 
- -bc_ - _comp Overrides the boundary components- - 
Calling Sequence of bcFunc and bcFunc_t#
If the type is DM_BC_ESSENTIAL
  void bcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[])
If the type is DM_BC_ESSENTIAL_FIELD or other _FIELD value,
  void bcFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
              PetscReal time, const PetscReal x[], PetscScalar bcval[])
- dim - the spatial dimension 
- Nf - the number of fields 
- uOff - the offset into u[] and u_t[] for each field 
- uOff_x - the offset into u_x[] for each field 
- u - each field evaluated at the current point 
- u_t - the time derivative of each field evaluated at the current point 
- u_x - the gradient of each field evaluated at the current point 
- aOff - the offset into a[] and a_t[] for each auxiliary field 
- aOff_x - the offset into a_x[] for each auxiliary field 
- a - each auxiliary field evaluated at the current point 
- a_t - the time derivative of each auxiliary field evaluated at the current point 
- a_x - the gradient of auxiliary each field evaluated at the current point 
- t - current time 
- x - coordinates of the current point 
- numConstants - number of constant parameters 
- constants - constant parameters 
- bcval - output values at the current point 
Notes#
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().
This function should only be used with DMFOREST currently, since labels cannot be defined before the underlying DMPLEX is built.
See Also#
PetscDS, PetscWeakForm, DMLabel, DMBoundaryConditionType, PetscDSAddBoundary(), PetscDSGetBoundary(), PetscDSSetResidual(), PetscDSSetBdResidual()
Level#
developer
Location#
Examples#
Index of all DT routines
Table of Contents for all manual pages
Index of all manual pages