DMAddBoundary#
Add a boundary condition to a model represented by a DM
Synopsis#
#include "petscdm.h"
#include "petscdmlabel.h"
#include "petscds.h"
PetscErrorCode DMAddBoundary(DM dm, 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, PetscInt *bd)
Collective
Input Parameters#
dm - The
DM
, with aPetscDS
that matches the problem being constrainedtype - The type of condition, e.g.
DM_BC_ESSENTIAL_ANALYTIC
,DM_BC_ESSENTIAL_FIELD
(Dirichlet), orDM_BC_NATURAL
(Neumann)name - The BC name
label - The label defining constrained points
Nv - The number of
DMLabel
values for constrained pointsvalues - An array of 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 deriative of the boundary values, or NULL
ctx - An optional user context for bcFunc
Output Parameter#
bd - (Optional) Boundary number
Options Database Keys#
-bc_
- -bc_
_comp Overrides the boundary components-
Notes#
Both bcFunc and bcFunc_t will depend on the boundary condition type. If the type if DM_BC_ESSENTIAL
, then the calling sequence is:
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, then the calling sequence is:
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
See Also#
DM Basics, DM
, DSGetBoundary()
, PetscDSAddBoundary()
Level#
intermediate
Location#
Examples#
src/tao/tutorials/ex2.c
src/snes/tutorials/ex62.c
src/snes/tutorials/ex20.c
src/snes/tutorials/ex27.c
src/tao/tutorials/ex1.c
src/snes/tutorials/ex71.c
src/snes/tutorials/ex63.c
src/snes/tutorials/ex56.c
src/snes/tutorials/ex17.c
src/snes/tutorials/ex23.c
Index of all DM routines
Table of Contents for all manual pages
Index of all manual pages