Actual source code: plexlabel.c

petsc-3.4.5 2014-06-29
  1: #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/

  5: PetscErrorCode DMLabelCreate(const char name[], DMLabel *label)
  6: {

 10:   PetscNew(struct _n_DMLabel, label);
 11:   PetscStrallocpy(name, &(*label)->name);

 13:   (*label)->refct          = 1;
 14:   (*label)->stratumValues  = NULL;
 15:   (*label)->stratumOffsets = NULL;
 16:   (*label)->stratumSizes   = NULL;
 17:   (*label)->points         = NULL;
 18:   (*label)->next           = NULL;
 19:   (*label)->pStart         = -1;
 20:   (*label)->pEnd           = -1;
 21:   (*label)->bt             = NULL;
 22:   return(0);
 23: }

 27: PetscErrorCode DMLabelGetName(DMLabel label, const char **name)
 28: {
 31:   *name = label->name;
 32:   return(0);
 33: }

 37: static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer)
 38: {
 39:   PetscInt       v;
 40:   PetscMPIInt    rank;

 44:   MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);
 45:   if (label) {
 46:     PetscViewerASCIIPrintf(viewer, "Label '%s':\n", label->name);
 47:     if (label->bt) {PetscViewerASCIIPrintf(viewer, "  Index has been calculated in [%d, %d)\n", label->pStart, label->pEnd);}
 48:     for (v = 0; v < label->numStrata; ++v) {
 49:       const PetscInt value = label->stratumValues[v];
 50:       PetscInt       p;

 52:       for (p = label->stratumOffsets[v]; p < label->stratumOffsets[v]+label->stratumSizes[v]; ++p) {
 53:         PetscViewerASCIISynchronizedPrintf(viewer, "[%D]: %D (%D)\n", rank, label->points[p], value);
 54:       }
 55:     }
 56:   }
 57:   PetscViewerFlush(viewer);
 58:   return(0);
 59: }

 63: PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
 64: {
 65:   PetscBool      iascii;

 70:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
 71:   if (iascii) {
 72:     DMLabelView_Ascii(label, viewer);
 73:   } else SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Viewer type %s not supported by this mesh object", ((PetscObject) viewer)->type_name);
 74:   return(0);
 75: }

 79: PetscErrorCode DMLabelDestroy(DMLabel *label)
 80: {

 84:   if (!(*label)) return(0);
 85:   if (--(*label)->refct > 0) return(0);
 86:   PetscFree((*label)->name);
 87:   PetscFree3((*label)->stratumValues,(*label)->stratumOffsets,(*label)->stratumSizes);
 88:   PetscFree((*label)->points);
 89:   PetscBTDestroy(&(*label)->bt);
 90:   PetscFree(*label);
 91:   return(0);
 92: }

 96: PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
 97: {
 98:   PetscInt       v, q;

102:   PetscNew(struct _n_DMLabel, labelnew);
103:   PetscStrallocpy(label->name, &(*labelnew)->name);

105:   (*labelnew)->refct     = 1;
106:   (*labelnew)->numStrata = label->numStrata;
107:   if (label->numStrata) {
108:     PetscMalloc3(label->numStrata,PetscInt,&(*labelnew)->stratumValues,label->numStrata+1,PetscInt,&(*labelnew)->stratumOffsets,label->numStrata,PetscInt,&(*labelnew)->stratumSizes);
109:     PetscMalloc(label->stratumOffsets[label->numStrata] * sizeof(PetscInt), &(*labelnew)->points);
110:     /* Could eliminate unused space here */
111:     for (v = 0; v < label->numStrata; ++v) {
112:       (*labelnew)->stratumValues[v]  = label->stratumValues[v];
113:       (*labelnew)->stratumOffsets[v] = label->stratumOffsets[v];
114:       (*labelnew)->stratumSizes[v]   = label->stratumSizes[v];
115:       for (q = label->stratumOffsets[v]; q < label->stratumOffsets[v]+label->stratumSizes[v]; ++q) {
116:         (*labelnew)->points[q] = label->points[q];
117:       }
118:     }
119:     (*labelnew)->stratumOffsets[label->numStrata] = label->stratumOffsets[label->numStrata];
120:   }
121:   (*labelnew)->pStart         = -1;
122:   (*labelnew)->pEnd           = -1;
123:   (*labelnew)->bt             = NULL;
124:   return(0);
125: }

129: /* This can be hooked into SetValue(),  ClearValue(), etc. for updating */
130: PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
131: {
132:   PetscInt       v;

136:   if (label->bt) {PetscBTDestroy(&label->bt);}
137:   label->pStart = pStart;
138:   label->pEnd   = pEnd;
139:   PetscBTCreate(pEnd - pStart, &label->bt);
140:   PetscBTMemzero(pEnd - pStart, label->bt);
141:   for (v = 0; v < label->numStrata; ++v) {
142:     PetscInt i;

144:     for (i = 0; i < label->stratumSizes[v]; ++i) {
145:       const PetscInt point = label->points[label->stratumOffsets[v]+i];

147:       if ((point < pStart) || (point >= pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %d is not in [%d, %d)", point, pStart, pEnd);
148:       PetscBTSet(label->bt, point);
149:     }
150:   }
151:   return(0);
152: }

156: PetscErrorCode DMLabelDestroyIndex(DMLabel label)
157: {

161:   label->pStart = -1;
162:   label->pEnd   = -1;
163:   if (label->bt) {PetscBTDestroy(&label->bt);}
164:   return(0);
165: }

169: /*@
170:   DMLabelHasPoint - Determine whether a label assigns a value to a point

172:   Input Parameters:
173: + label - the DMLabel
174: - point - the point

176:   Output Parameter:
177: . contains - Flag indicating whether the label maps this point to a value

179:   Note: The user must call DMLabelCreateIndex() before this function.

181:   Level: developer

183: .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
184: @*/
185: PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
186: {
189: #if defined(PETSC_USE_DEBUG)
190:   if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
191:   if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %d is not in [%d, %d)", point, label->pStart, label->pEnd);
192: #endif
193:   *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
194:   return(0);
195: }

199: PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
200: {
201:   PetscInt       v;

206:   *value = -1;
207:   for (v = 0; v < label->numStrata; ++v) {
208:     PetscInt i;

210:     PetscFindInt(point, label->stratumSizes[v], &label->points[label->stratumOffsets[v]], &i);
211:     if (i >= 0) {
212:       *value = label->stratumValues[v];
213:       break;
214:     }
215:   }
216:   return(0);
217: }

221: PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
222: {
223:   PetscInt       v, loc;

227:   /* Find, or add, label value */
228:   for (v = 0; v < label->numStrata; ++v) {
229:     if (label->stratumValues[v] == value) break;
230:   }
231:   if (v >= label->numStrata) {
232:     PetscInt *tmpV, *tmpO, *tmpS;

234:     PetscMalloc3(label->numStrata+1,PetscInt,&tmpV,label->numStrata+2,PetscInt,&tmpO,label->numStrata+1,PetscInt,&tmpS);
235:     for (v = 0; v < label->numStrata; ++v) {
236:       tmpV[v] = label->stratumValues[v];
237:       tmpO[v] = label->stratumOffsets[v];
238:       tmpS[v] = label->stratumSizes[v];
239:     }
240:     tmpV[v]   = value;
241:     tmpO[v]   = v == 0 ? 0 : label->stratumOffsets[v];
242:     tmpS[v]   = 0;
243:     tmpO[v+1] = tmpO[v];
244:     ++label->numStrata;
245:     PetscFree3(label->stratumValues,label->stratumOffsets,label->stratumSizes);

247:     label->stratumValues  = tmpV;
248:     label->stratumOffsets = tmpO;
249:     label->stratumSizes   = tmpS;
250:   }
251:   /* Check whether point exists */
252:   PetscFindInt(point,label->stratumSizes[v],label->points+label->stratumOffsets[v],&loc);
253:   if (loc < 0) {
254:     PetscInt off = label->stratumOffsets[v] - (loc+1); /* decode insert location */
255:     /* Check for reallocation */
256:     if (label->stratumSizes[v] >= label->stratumOffsets[v+1]-label->stratumOffsets[v]) {
257:       PetscInt oldSize   = label->stratumOffsets[v+1]-label->stratumOffsets[v];
258:       PetscInt newSize   = PetscMax(10, 2*oldSize);  /* Double the size, since 2 is the optimal base for this online algorithm */
259:       PetscInt shift     = newSize - oldSize;
260:       PetscInt allocSize = label->stratumOffsets[label->numStrata] + shift;
261:       PetscInt *newPoints;
262:       PetscInt w, q;

264:       PetscMalloc(allocSize * sizeof(PetscInt), &newPoints);
265:       for (q = 0; q < label->stratumOffsets[v]+label->stratumSizes[v]; ++q) {
266:         newPoints[q] = label->points[q];
267:       }
268:       for (w = v+1; w < label->numStrata; ++w) {
269:         for (q = label->stratumOffsets[w]; q < label->stratumOffsets[w]+label->stratumSizes[w]; ++q) {
270:           newPoints[q+shift] = label->points[q];
271:         }
272:         label->stratumOffsets[w] += shift;
273:       }
274:       label->stratumOffsets[label->numStrata] += shift;

276:       PetscFree(label->points);
277:       label->points = newPoints;
278:     }
279:     PetscMemmove(&label->points[off+1], &label->points[off], (label->stratumSizes[v]+(loc+1)) * sizeof(PetscInt));

281:     label->points[off] = point;
282:     ++label->stratumSizes[v];
283:     if (label->bt) {
284:       if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %d is not in [%d, %d)", point, label->pStart, label->pEnd);
285:       PetscBTSet(label->bt, point);
286:     }
287:   }
288:   return(0);
289: }

293: PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
294: {
295:   PetscInt       v, p;

299:   /* Find label value */
300:   for (v = 0; v < label->numStrata; ++v) {
301:     if (label->stratumValues[v] == value) break;
302:   }
303:   if (v >= label->numStrata) return(0);
304:   /* Check whether point exists */
305:   for (p = label->stratumOffsets[v]; p < label->stratumOffsets[v]+label->stratumSizes[v]; ++p) {
306:     if (label->points[p] == point) {
307:       /* Found point */
308:       PetscInt q;

310:       for (q = p+1; q < label->stratumOffsets[v]+label->stratumSizes[v]; ++q) {
311:         label->points[q-1] = label->points[q];
312:       }
313:       --label->stratumSizes[v];
314:       if (label->bt) {
315:         if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %d is not in [%d, %d)", point, label->pStart, label->pEnd);
316:         PetscBTClear(label->bt, point);
317:       }
318:       break;
319:     }
320:   }
321:   return(0);
322: }

326: PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
327: {
330:   *numValues = label->numStrata;
331:   return(0);
332: }

336: PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
337: {

342:   ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_COPY_VALUES, values);
343:   return(0);
344: }

348: PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
349: {
350:   PetscInt v;

354:   *size = 0;
355:   for (v = 0; v < label->numStrata; ++v) {
356:     if (label->stratumValues[v] == value) {
357:       *size = label->stratumSizes[v];
358:       break;
359:     }
360:   }
361:   return(0);
362: }

366: PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
367: {
368:   PetscInt       v;

373:   *points = NULL;
374:   for (v = 0; v < label->numStrata; ++v) {
375:     if (label->stratumValues[v] == value) {
376:       ISCreateGeneral(PETSC_COMM_SELF, label->stratumSizes[v], &label->points[label->stratumOffsets[v]], PETSC_COPY_VALUES, points);
377:       break;
378:     }
379:   }
380:   return(0);
381: }

385: PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
386: {
387:   PetscInt       v;

391:   for (v = 0; v < label->numStrata; ++v) {
392:     if (label->stratumValues[v] == value) break;
393:   }
394:   if (v >= label->numStrata) return(0);
395:   if (label->bt) {
396:     PetscInt i;

398:     for (i = 0; i < label->stratumSizes[v]; ++i) {
399:       const PetscInt point = label->points[label->stratumOffsets[v]+i];

401:       if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %d is not in [%d, %d)", point, label->pStart, label->pEnd);
402:       PetscBTClear(label->bt, point);
403:     }
404:   }
405:   label->stratumSizes[v] = 0;
406:   return(0);
407: }



413: /*@C
414:   DMPlexCreateLabel - Create a label of the given name if it does not already exist

416:   Not Collective

418:   Input Parameters:
419: + dm   - The DMPlex object
420: - name - The label name

422:   Level: intermediate

424: .keywords: mesh
425: .seealso: DMLabelCreate(), DMPlexHasLabel(), DMPlexGetLabelValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
426: @*/
427: PetscErrorCode DMPlexCreateLabel(DM dm, const char name[])
428: {
429:   DM_Plex        *mesh = (DM_Plex*) dm->data;
430:   DMLabel        next  = mesh->labels;
431:   PetscBool      flg   = PETSC_FALSE;

437:   while (next) {
438:     PetscStrcmp(name, next->name, &flg);
439:     if (flg) break;
440:     next = next->next;
441:   }
442:   if (!flg) {
443:     DMLabel tmpLabel = mesh->labels;

445:     DMLabelCreate(name, &mesh->labels);

447:     mesh->labels->next = tmpLabel;
448:   }
449:   return(0);
450: }

454: /*@C
455:   DMPlexGetLabelValue - Get the value in a Sieve Label for the given point, with 0 as the default

457:   Not Collective

459:   Input Parameters:
460: + dm   - The DMPlex object
461: . name - The label name
462: - point - The mesh point

464:   Output Parameter:
465: . value - The label value for this point, or -1 if the point is not in the label

467:   Level: beginner

469: .keywords: mesh
470: .seealso: DMLabelGetValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
471: @*/
472: PetscErrorCode DMPlexGetLabelValue(DM dm, const char name[], PetscInt point, PetscInt *value)
473: {
474:   DMLabel        label;

480:   DMPlexGetLabel(dm, name, &label);
481:   if (!label) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No label named %s was found", name);
482:   DMLabelGetValue(label, point, value);
483:   return(0);
484: }

488: /*@C
489:   DMPlexSetLabelValue - Add a point to a Sieve Label with given value

491:   Not Collective

493:   Input Parameters:
494: + dm   - The DMPlex object
495: . name - The label name
496: . point - The mesh point
497: - value - The label value for this point

499:   Output Parameter:

501:   Level: beginner

503: .keywords: mesh
504: .seealso: DMLabelSetValue(), DMPlexGetStratumIS(), DMPlexClearLabelValue()
505: @*/
506: PetscErrorCode DMPlexSetLabelValue(DM dm, const char name[], PetscInt point, PetscInt value)
507: {
508:   DMLabel        label;

514:   DMPlexGetLabel(dm, name, &label);
515:   if (!label) {
516:     DMPlexCreateLabel(dm, name);
517:     DMPlexGetLabel(dm, name, &label);
518:   }
519:   DMLabelSetValue(label, point, value);
520:   return(0);
521: }

525: /*@C
526:   DMPlexClearLabelValue - Remove a point from a Sieve Label with given value

528:   Not Collective

530:   Input Parameters:
531: + dm   - The DMPlex object
532: . name - The label name
533: . point - The mesh point
534: - value - The label value for this point

536:   Output Parameter:

538:   Level: beginner

540: .keywords: mesh
541: .seealso: DMLabelClearValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
542: @*/
543: PetscErrorCode DMPlexClearLabelValue(DM dm, const char name[], PetscInt point, PetscInt value)
544: {
545:   DMLabel        label;

551:   DMPlexGetLabel(dm, name, &label);
552:   if (!label) return(0);
553:   DMLabelClearValue(label, point, value);
554:   return(0);
555: }

559: /*@C
560:   DMPlexGetLabelSize - Get the number of different integer ids in a Label

562:   Not Collective

564:   Input Parameters:
565: + dm   - The DMPlex object
566: - name - The label name

568:   Output Parameter:
569: . size - The number of different integer ids, or 0 if the label does not exist

571:   Level: beginner

573: .keywords: mesh
574: .seealso: DMLabeGetNumValues(), DMPlexSetLabelValue()
575: @*/
576: PetscErrorCode DMPlexGetLabelSize(DM dm, const char name[], PetscInt *size)
577: {
578:   DMLabel        label;

585:   DMPlexGetLabel(dm, name, &label);
586:   *size = 0;
587:   if (!label) return(0);
588:   DMLabelGetNumValues(label, size);
589:   return(0);
590: }

594: /*@C
595:   DMPlexGetLabelIdIS - Get the integer ids in a label

597:   Not Collective

599:   Input Parameters:
600: + mesh - The DMPlex object
601: - name - The label name

603:   Output Parameter:
604: . ids - The integer ids, or NULL if the label does not exist

606:   Level: beginner

608: .keywords: mesh
609: .seealso: DMLabelGetValueIS(), DMPlexGetLabelSize()
610: @*/
611: PetscErrorCode DMPlexGetLabelIdIS(DM dm, const char name[], IS *ids)
612: {
613:   DMLabel        label;

620:   DMPlexGetLabel(dm, name, &label);
621:   *ids = NULL;
622:   if (!label) return(0);
623:   DMLabelGetValueIS(label, ids);
624:   return(0);
625: }

629: /*@C
630:   DMPlexGetStratumSize - Get the number of points in a label stratum

632:   Not Collective

634:   Input Parameters:
635: + dm - The DMPlex object
636: . name - The label name
637: - value - The stratum value

639:   Output Parameter:
640: . size - The stratum size

642:   Level: beginner

644: .keywords: mesh
645: .seealso: DMLabelGetStratumSize(), DMPlexGetLabelSize(), DMPlexGetLabelIds()
646: @*/
647: PetscErrorCode DMPlexGetStratumSize(DM dm, const char name[], PetscInt value, PetscInt *size)
648: {
649:   DMLabel        label;

656:   DMPlexGetLabel(dm, name, &label);
657:   *size = 0;
658:   if (!label) return(0);
659:   DMLabelGetStratumSize(label, value, size);
660:   return(0);
661: }

665: /*@C
666:   DMPlexGetStratumIS - Get the points in a label stratum

668:   Not Collective

670:   Input Parameters:
671: + dm - The DMPlex object
672: . name - The label name
673: - value - The stratum value

675:   Output Parameter:
676: . points - The stratum points, or NULL if the label does not exist or does not have that value

678:   Level: beginner

680: .keywords: mesh
681: .seealso: DMLabelGetStratumIS(), DMPlexGetStratumSize()
682: @*/
683: PetscErrorCode DMPlexGetStratumIS(DM dm, const char name[], PetscInt value, IS *points)
684: {
685:   DMLabel        label;

692:   DMPlexGetLabel(dm, name, &label);
693:   *points = NULL;
694:   if (!label) return(0);
695:   DMLabelGetStratumIS(label, value, points);
696:   return(0);
697: }

701: /*@C
702:   DMPlexClearLabelStratum - Remove all points from a stratum from a Sieve Label

704:   Not Collective

706:   Input Parameters:
707: + dm   - The DMPlex object
708: . name - The label name
709: - value - The label value for this point

711:   Output Parameter:

713:   Level: beginner

715: .keywords: mesh
716: .seealso: DMLabelClearStratum(), DMPlexSetLabelValue(), DMPlexGetStratumIS(), DMPlexClearLabelValue()
717: @*/
718: PetscErrorCode DMPlexClearLabelStratum(DM dm, const char name[], PetscInt value)
719: {
720:   DMLabel        label;

726:   DMPlexGetLabel(dm, name, &label);
727:   if (!label) return(0);
728:   DMLabelClearStratum(label, value);
729:   return(0);
730: }



736: /*@
737:   DMPlexGetNumLabels - Return the number of labels defined by the mesh

739:   Not Collective

741:   Input Parameter:
742: . dm   - The DMPlex object

744:   Output Parameter:
745: . numLabels - the number of Labels

747:   Level: intermediate

749: .keywords: mesh
750: .seealso: DMPlexGetLabelValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
751: @*/
752: PetscErrorCode DMPlexGetNumLabels(DM dm, PetscInt *numLabels)
753: {
754:   DM_Plex  *mesh = (DM_Plex*) dm->data;
755:   DMLabel  next  = mesh->labels;
756:   PetscInt n     = 0;

761:   while (next) {
762:     ++n;
763:     next = next->next;
764:   }
765:   *numLabels = n;
766:   return(0);
767: }

771: /*@C
772:   DMPlexGetLabelName - Return the name of nth label

774:   Not Collective

776:   Input Parameters:
777: + dm - The DMPlex object
778: - n  - the label number

780:   Output Parameter:
781: . name - the label name

783:   Level: intermediate

785: .keywords: mesh
786: .seealso: DMPlexGetLabelValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
787: @*/
788: PetscErrorCode DMPlexGetLabelName(DM dm, PetscInt n, const char **name)
789: {
790:   DM_Plex  *mesh = (DM_Plex*) dm->data;
791:   DMLabel  next  = mesh->labels;
792:   PetscInt l     = 0;

797:   while (next) {
798:     if (l == n) {
799:       *name = next->name;
800:       return(0);
801:     }
802:     ++l;
803:     next = next->next;
804:   }
805:   SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %d does not exist in this DM", n);
806: }

810: /*@C
811:   DMPlexHasLabel - Determine whether the mesh has a label of a given name

813:   Not Collective

815:   Input Parameters:
816: + dm   - The DMPlex object
817: - name - The label name

819:   Output Parameter:
820: . hasLabel - PETSC_TRUE if the label is present

822:   Level: intermediate

824: .keywords: mesh
825: .seealso: DMPlexCreateLabel(), DMPlexGetLabelValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
826: @*/
827: PetscErrorCode DMPlexHasLabel(DM dm, const char name[], PetscBool *hasLabel)
828: {
829:   DM_Plex        *mesh = (DM_Plex*) dm->data;
830:   DMLabel        next  = mesh->labels;

837:   *hasLabel = PETSC_FALSE;
838:   while (next) {
839:     PetscStrcmp(name, next->name, hasLabel);
840:     if (*hasLabel) break;
841:     next = next->next;
842:   }
843:   return(0);
844: }

848: /*@C
849:   DMPlexGetLabel - Return the label of a given name, or NULL

851:   Not Collective

853:   Input Parameters:
854: + dm   - The DMPlex object
855: - name - The label name

857:   Output Parameter:
858: . label - The DMLabel, or NULL if the label is absent

860:   Level: intermediate

862: .keywords: mesh
863: .seealso: DMPlexCreateLabel(), DMPlexHasLabel(), DMPlexGetLabelValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
864: @*/
865: PetscErrorCode DMPlexGetLabel(DM dm, const char name[], DMLabel *label)
866: {
867:   DM_Plex        *mesh = (DM_Plex*) dm->data;
868:   DMLabel        next  = mesh->labels;
869:   PetscBool      hasLabel;

876:   *label = NULL;
877:   while (next) {
878:     PetscStrcmp(name, next->name, &hasLabel);
879:     if (hasLabel) {
880:       *label = next;
881:       break;
882:     }
883:     next = next->next;
884:   }
885:   return(0);
886: }

890: /*@C
891:   DMPlexAddLabel - Add the label to this mesh

893:   Not Collective

895:   Input Parameters:
896: + dm   - The DMPlex object
897: - label - The DMLabel

899:   Level: developer

901: .keywords: mesh
902: .seealso: DMPlexCreateLabel(), DMPlexHasLabel(), DMPlexGetLabelValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
903: @*/
904: PetscErrorCode DMPlexAddLabel(DM dm, DMLabel label)
905: {
906:   DM_Plex        *mesh = (DM_Plex*) dm->data;
907:   PetscBool      hasLabel;

912:   DMPlexHasLabel(dm, label->name, &hasLabel);
913:   if (hasLabel) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %s already exists in this DM", label->name);
914:   label->next  = mesh->labels;
915:   mesh->labels = label;
916:   return(0);
917: }

921: /*@C
922:   DMPlexRemoveLabel - Remove the label from this mesh

924:   Not Collective

926:   Input Parameters:
927: + dm   - The DMPlex object
928: - name - The label name

930:   Output Parameter:
931: . label - The DMLabel, or NULL if the label is absent

933:   Level: developer

935: .keywords: mesh
936: .seealso: DMPlexCreateLabel(), DMPlexHasLabel(), DMPlexGetLabelValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
937: @*/
938: PetscErrorCode DMPlexRemoveLabel(DM dm, const char name[], DMLabel *label)
939: {
940:   DM_Plex        *mesh = (DM_Plex*) dm->data;
941:   DMLabel        next  = mesh->labels;
942:   DMLabel        last  = NULL;
943:   PetscBool      hasLabel;

948:   DMPlexHasLabel(dm, name, &hasLabel);
949:   *label = NULL;
950:   if (!hasLabel) return(0);
951:   while (next) {
952:     PetscStrcmp(name, next->name, &hasLabel);
953:     if (hasLabel) {
954:       if (last) last->next   = next->next;
955:       else      mesh->labels = next->next;
956:       next->next = NULL;
957:       *label     = next;
958:       break;
959:     }
960:     last = next;
961:     next = next->next;
962:   }
963:   return(0);
964: }