Actual source code: plexbc.c
petsc-3.6.1 2015-08-06
1: #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/
5: static PetscErrorCode BoundaryDuplicate(DMBoundary bd, DMBoundary *boundary)
6: {
7: DMBoundary b = bd, b2, bold = NULL;
11: *boundary = NULL;
12: for (; b; b = b->next, bold = b2) {
13: PetscNew(&b2);
14: PetscStrallocpy(b->name, (char **) &b2->name);
15: PetscStrallocpy(b->labelname, (char **) &b2->labelname);
16: PetscMalloc1(b->numids, &b2->ids);
17: PetscMemcpy(b2->ids, b->ids, b->numids*sizeof(PetscInt));
18: PetscMalloc1(b->numcomps, &b2->comps);
19: PetscMemcpy(b2->comps, b->comps, b->numcomps*sizeof(PetscInt));
20: b2->label = NULL;
21: b2->essential = b->essential;
22: b2->field = b->field;
23: b2->numcomps = b->numcomps;
24: b2->func = b->func;
25: b2->numids = b->numids;
26: b2->ctx = b->ctx;
27: b2->next = NULL;
28: if (!*boundary) *boundary = b2;
29: if (bold) bold->next = b2;
30: }
31: return(0);
32: }
36: PetscErrorCode DMPlexCopyBoundary(DM dm, DM dmNew)
37: {
38: DM_Plex *mesh = (DM_Plex *) dm->data;
39: DM_Plex *meshNew = (DM_Plex *) dmNew->data;
40: DMBoundary b;
44: BoundaryDuplicate(mesh->boundary, &meshNew->boundary);
45: for (b = meshNew->boundary; b; b = b->next) {
46: if (b->labelname) {
47: DMPlexGetLabel(dmNew, b->labelname, &b->label);
48: if (!b->label) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Label %s does not exist in this DM", b->labelname);
49: }
50: }
51: return(0);
52: }
56: /*@C
57: DMPlexAddBoundary - Add a boundary condition to the model
59: Input Parameters:
60: + dm - The mesh object
61: . isEssential - Flag for an essential (Dirichlet) condition, as opposed to a natural (Neumann) condition
62: . name - The BC name
63: . labelname - The label defining constrained points
64: . field - The field to constrain
65: . numcomps - The number of constrained field components
66: . comps - An array of constrained component numbers
67: . bcFunc - A pointwise function giving boundary values
68: . numids - The number of DMLabel ids for constrained points
69: . ids - An array of ids for constrained points
70: - ctx - An optional user context for bcFunc
72: Options Database Keys:
73: + -bc_<boundary name> <num> - Overrides the boundary ids
74: - -bc_<boundary name>_comp <num> - Overrides the boundary components
76: Level: developer
78: .seealso: DMPlexGetBoundary()
79: @*/
80: PetscErrorCode DMPlexAddBoundary(DM dm, PetscBool isEssential, const char name[], const char labelname[], PetscInt field, PetscInt numcomps, const PetscInt *comps, void (*bcFunc)(), PetscInt numids, const PetscInt *ids, void *ctx)
81: {
82: DM_Plex *mesh = (DM_Plex *) dm->data;
83: DMBoundary b;
88: PetscNew(&b);
89: PetscStrallocpy(name, (char **) &b->name);
90: PetscStrallocpy(labelname, (char **) &b->labelname);
91: PetscMalloc1(numcomps, &b->comps);
92: if (numcomps) {PetscMemcpy(b->comps, comps, numcomps*sizeof(PetscInt));}
93: PetscMalloc1(numids, &b->ids);
94: if (numids) {PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));}
95: if (b->labelname) {
96: DMPlexGetLabel(dm, b->labelname, &b->label);
97: if (!b->label) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Label %s does not exist in this DM", b->labelname);
98: }
99: b->essential = isEssential;
100: b->field = field;
101: b->numcomps = numcomps;
102: b->func = bcFunc;
103: b->numids = numids;
104: b->ctx = ctx;
105: b->next = mesh->boundary;
106: mesh->boundary = b;
107: return(0);
108: }
112: /*@
113: DMPlexGetNumBoundary - Get the number of registered BC
115: Input Parameters:
116: . dm - The mesh object
118: Output Parameters:
119: . numBd - The number of BC
121: Level: intermediate
123: .seealso: DMPlexAddBoundary(), DMPlexGetBoundary()
124: @*/
125: PetscErrorCode DMPlexGetNumBoundary(DM dm, PetscInt *numBd)
126: {
127: DM_Plex *mesh = (DM_Plex *) dm->data;
128: DMBoundary b = mesh->boundary;
133: *numBd = 0;
134: while (b) {++(*numBd); b = b->next;}
135: return(0);
136: }
140: /*@C
141: DMPlexAddBoundary - Add a boundary condition to the model
143: Input Parameters:
144: + dm - The mesh object
145: - bd - The BC number
147: Output Parameters:
148: + isEssential - Flag for an essential (Dirichlet) condition, as opposed to a natural (Neumann) condition
149: . name - The BC name
150: . labelname - The label defining constrained points
151: . field - The field to constrain
152: . numcomps - The number of constrained field components
153: . comps - An array of constrained component numbers
154: . bcFunc - A pointwise function giving boundary values
155: . numids - The number of DMLabel ids for constrained points
156: . ids - An array of ids for constrained points
157: - ctx - An optional user context for bcFunc
159: Options Database Keys:
160: + -bc_<boundary name> <num> - Overrides the boundary ids
161: - -bc_<boundary name>_comp <num> - Overrides the boundary components
163: Level: developer
165: .seealso: DMPlexAddBoundary()
166: @*/
167: PetscErrorCode DMPlexGetBoundary(DM dm, PetscInt bd, PetscBool *isEssential, const char **name, const char **labelname, PetscInt *field, PetscInt *numcomps, const PetscInt **comps, void (**func)(), PetscInt *numids, const PetscInt **ids, void **ctx)
168: {
169: DM_Plex *mesh = (DM_Plex *) dm->data;
170: DMBoundary b = mesh->boundary;
171: PetscInt n = 0;
175: while (b) {
176: if (n == bd) break;
177: b = b->next;
178: ++n;
179: }
180: if (n != bd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n);
181: if (isEssential) {
183: *isEssential = b->essential;
184: }
185: if (name) {
187: *name = b->name;
188: }
189: if (labelname) {
191: *labelname = b->labelname;
192: }
193: if (field) {
195: *field = b->field;
196: }
197: if (numcomps) {
199: *numcomps = b->numcomps;
200: }
201: if (comps) {
203: *comps = b->comps;
204: }
205: if (func) {
207: *func = b->func;
208: }
209: if (numids) {
211: *numids = b->numids;
212: }
213: if (ids) {
215: *ids = b->ids;
216: }
217: if (ctx) {
219: *ctx = b->ctx;
220: }
221: return(0);
222: }
226: PetscErrorCode DMPlexIsBoundaryPoint(DM dm, PetscInt point, PetscBool *isBd)
227: {
228: DM_Plex *mesh = (DM_Plex *) dm->data;
229: DMBoundary b = mesh->boundary;
235: *isBd = PETSC_FALSE;
236: while (b && !(*isBd)) {
237: if (b->label) {
238: PetscInt i;
240: for (i = 0; i < b->numids && !(*isBd); ++i) {
241: DMLabelStratumHasPoint(b->label, b->ids[i], point, isBd);
242: }
243: }
244: b = b->next;
245: }
246: return(0);
247: }