petsc-3.14.5 2021-03-03
DMAddBoundary
Add a boundary condition to the model
Synopsis
#include "petscdm.h"
#include "petscdmlabel.h"
#include "petscds.h"
PetscErrorCode DMAddBoundary(DM dm, DMBoundaryConditionType type, const char name[], const char labelname[], PetscInt field, PetscInt numcomps, const PetscInt *comps, void (*bcFunc)(void), void (*bcFunc_t)(void), PetscInt numids, const PetscInt *ids, void *ctx)
Collective on dm
Input Parameters
| dm | - The DM, with a PetscDS that matches the problem being constrained
|
| type | - The type of condition, e.g. DM_BC_ESSENTIAL_ANALYTIC/DM_BC_ESSENTIAL_FIELD (Dirichlet), or DM_BC_NATURAL (Neumann)
|
| name | - The BC name
|
| labelname | - The label defining constrained points
|
| field | - The field to constrain
|
| numcomps | - 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
|
| numids | - The number of DMLabel ids for constrained points
|
| ids | - An array of ids for constrained points
|
| ctx | - An optional user context for bcFunc
|
Options Database Keys
| -bc_<boundary name> <num> | - Overrides the boundary ids
|
| -bc_<boundary name>_comp <num> | - Overrides the boundary components
|
Note
Both bcFunc abd bcFunc_t will depend on the boundary condition type. If the type if DM_BC_ESSENTIAL, Then the calling sequence is
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
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
DMGetBoundary(), PetscDSAddBoundary()
Level
developer
Location
src/dm/interface/dm.c
Examples
src/dm/impls/plex/tutorials/ex2.c.html
src/snes/tutorials/ex12.c.html
src/snes/tutorials/ex17.c.html
src/snes/tutorials/ex56.c.html
src/snes/tutorials/ex62.c.html
src/snes/tutorials/ex77.c.html
src/ts/tutorials/ex45.c.html
src/ts/tutorials/ex46.c.html
src/ts/tutorials/ex48.c.html
src/tao/tutorials/ex1.c.html
src/tao/tutorials/ex2.c.html
Index of all DM routines
Table of Contents for all manual pages
Index of all manual pages