Actual source code: plexlabel.c

petsc-dev 2014-02-02
Report Typos and Errors
  1: #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/

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

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

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

 30: static PetscErrorCode DMLabelMakeValid_Private(DMLabel label)
 31: {
 32:   PetscInt       off = 0, v;

 35:   if (label->arrayValid) return 0;
 37:   PetscMalloc2(label->numStrata,&label->stratumSizes,label->numStrata+1,&label->stratumOffsets);
 38:   for (v = 0; v < label->numStrata; ++v) {
 39:     PetscInt size = 0;
 40:     PetscHashISize(label->ht[v], size);
 41:     label->stratumSizes[v]   = size;
 42:     label->stratumOffsets[v] = off;
 43:     off += size;
 44:   }
 45:   label->stratumOffsets[v] = off;
 46:   PetscMalloc1(off, &label->points);
 47:   off = 0;
 48:   for (v = 0; v < label->numStrata; ++v) {
 49:     PetscHashIGetKeys(label->ht[v], &off, label->points);
 50:     if (off != label->stratumOffsets[v+1]) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of contributed points %d from value %d should be %d", off-label->stratumOffsets[v], label->stratumValues[v], label->stratumOffsets[v+1]-label->stratumOffsets[v]);
 51:     PetscHashIDestroy(label->ht[v]);
 52:     PetscSortInt(label->stratumSizes[v], &label->points[label->stratumOffsets[v]]);
 53:     if (label->bt) {
 54:       PetscInt p;

 56:       for (p = label->stratumOffsets[v]; p < label->stratumOffsets[v]+label->stratumSizes[v]; ++p) {
 57:         const PetscInt point = label->points[p];

 59:         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);
 60:         PetscBTSet(label->bt, point - label->pStart);
 61:       }
 62:     }
 63:   }
 64:   PetscFree(label->ht);
 65:   label->ht = NULL;
 66:   label->arrayValid = PETSC_TRUE;
 67:   return(0);
 68: }

 72: static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label)
 73: {
 74:   PetscInt       v;

 77:   if (!label->arrayValid) return 0;
 79:   PetscMalloc1(label->numStrata, &label->ht);
 80:   for (v = 0; v < label->numStrata; ++v) {
 81:     PETSC_UNUSED PetscHashIIter ret, iter;
 82:     PetscInt                    p;

 84:     PetscHashICreate(label->ht[v]);
 85:     for (p = label->stratumOffsets[v]; p < label->stratumOffsets[v]+label->stratumSizes[v]; ++p) PetscHashIPut(label->ht[v], label->points[p], ret, iter);
 86:   }
 87:   PetscFree2(label->stratumSizes,label->stratumOffsets);
 88:   PetscFree(label->points);
 89:   label->arrayValid = PETSC_FALSE;
 90:   return(0);
 91: }

 95: PetscErrorCode DMLabelGetName(DMLabel label, const char **name)
 96: {
 99:   *name = label->name;
100:   return(0);
101: }

105: static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer)
106: {
107:   PetscInt       v;
108:   PetscMPIInt    rank;

112:   MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);
113:   if (label) {
114:     PetscViewerASCIIPrintf(viewer, "Label '%s':\n", label->name);
115:     if (label->bt) {PetscViewerASCIIPrintf(viewer, "  Index has been calculated in [%d, %d)\n", label->pStart, label->pEnd);}
116:     for (v = 0; v < label->numStrata; ++v) {
117:       const PetscInt value = label->stratumValues[v];
118:       PetscInt       p;

120:       for (p = label->stratumOffsets[v]; p < label->stratumOffsets[v]+label->stratumSizes[v]; ++p) {
121:         PetscViewerASCIISynchronizedPrintf(viewer, "[%D]: %D (%D)\n", rank, label->points[p], value);
122:       }
123:     }
124:   }
125:   PetscViewerFlush(viewer);
126:   return(0);
127: }

131: PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
132: {
133:   PetscBool      iascii;

138:   if (label) {DMLabelMakeValid_Private(label);}
139:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
140:   if (iascii) {
141:     DMLabelView_Ascii(label, viewer);
142:   } else SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Viewer type %s not supported by this mesh object", ((PetscObject) viewer)->type_name);
143:   return(0);
144: }

148: PetscErrorCode DMLabelDestroy(DMLabel *label)
149: {

153:   if (!(*label)) return(0);
154:   if (--(*label)->refct > 0) return(0);
155:   PetscFree((*label)->name);
156:   PetscFree((*label)->stratumValues);
157:   PetscFree2((*label)->stratumSizes,(*label)->stratumOffsets);
158:   PetscFree((*label)->points);
159:   if ((*label)->ht) {
160:     PetscInt v;
161:     for (v = 0; v < (*label)->numStrata; ++v) {PetscHashIDestroy((*label)->ht[v]);}
162:     PetscFree((*label)->ht);
163:   }
164:   PetscBTDestroy(&(*label)->bt);
165:   PetscFree(*label);
166:   return(0);
167: }

171: PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
172: {
173:   PetscInt       v, q;

177:   DMLabelMakeValid_Private(label);
178:   PetscNew(labelnew);
179:   PetscStrallocpy(label->name, &(*labelnew)->name);

181:   (*labelnew)->refct      = 1;
182:   (*labelnew)->numStrata  = label->numStrata;
183:   (*labelnew)->arrayValid = PETSC_TRUE;
184:   if (label->numStrata) {
185:     PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);
186:     PetscMalloc2(label->numStrata,&(*labelnew)->stratumSizes,label->numStrata+1,&(*labelnew)->stratumOffsets);
187:     PetscMalloc1(label->stratumOffsets[label->numStrata], &(*labelnew)->points);
188:     /* Could eliminate unused space here */
189:     for (v = 0; v < label->numStrata; ++v) {
190:       (*labelnew)->stratumValues[v]  = label->stratumValues[v];
191:       (*labelnew)->stratumOffsets[v] = label->stratumOffsets[v];
192:       (*labelnew)->stratumSizes[v]   = label->stratumSizes[v];
193:       for (q = label->stratumOffsets[v]; q < label->stratumOffsets[v]+label->stratumSizes[v]; ++q) {
194:         (*labelnew)->points[q] = label->points[q];
195:       }
196:     }
197:     (*labelnew)->stratumOffsets[label->numStrata] = label->stratumOffsets[label->numStrata];
198:   }
199:   (*labelnew)->ht     = NULL;
200:   (*labelnew)->pStart = -1;
201:   (*labelnew)->pEnd   = -1;
202:   (*labelnew)->bt     = NULL;
203:   return(0);
204: }

208: /* This can be hooked into SetValue(),  ClearValue(), etc. for updating */
209: PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
210: {
211:   PetscInt       v;

215:   DMLabelMakeValid_Private(label);
216:   if (label->bt) {PetscBTDestroy(&label->bt);}
217:   label->pStart = pStart;
218:   label->pEnd   = pEnd;
219:   PetscBTCreate(pEnd - pStart, &label->bt);
220:   PetscBTMemzero(pEnd - pStart, label->bt);
221:   for (v = 0; v < label->numStrata; ++v) {
222:     PetscInt i;

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

227:       if ((point < pStart) || (point >= pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %d is not in [%d, %d)", point, pStart, pEnd);
228:       PetscBTSet(label->bt, point - pStart);
229:     }
230:   }
231:   return(0);
232: }

236: PetscErrorCode DMLabelDestroyIndex(DMLabel label)
237: {

241:   label->pStart = -1;
242:   label->pEnd   = -1;
243:   if (label->bt) {PetscBTDestroy(&label->bt);}
244:   return(0);
245: }

249: /*@
250:   DMLabelHasPoint - Determine whether a label assigns a value to a point

252:   Input Parameters:
253: + label - the DMLabel
254: - point - the point

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

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

261:   Level: developer

263: .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
264: @*/
265: PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
266: {

271:   DMLabelMakeValid_Private(label);
272: #if defined(PETSC_USE_DEBUG)
273:   if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
274:   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);
275: #endif
276:   *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
277:   return(0);
278: }

282: /*@
283:   DMLabelGetValue - Return the value a label assigns to a point, or -1

285:   Input Parameters:
286: + label - the DMLabel
287: - point - the point

289:   Output Parameter:
290: . value - The point value, or -1

292:   Level: intermediate

294: .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue()
295: @*/
296: PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
297: {
298:   PetscInt       v;

303:   *value = -1;
304:   for (v = 0; v < label->numStrata; ++v) {
305:     if (label->arrayValid) {
306:       PetscInt i;

308:       PetscFindInt(point, label->stratumSizes[v], &label->points[label->stratumOffsets[v]], &i);
309:       if (i >= 0) {
310:         *value = label->stratumValues[v];
311:         break;
312:       }
313:     } else {
314:       PetscBool has;

316:       PetscHashIHasKey(label->ht[v], point, has);
317:       if (has) {
318:         *value = label->stratumValues[v];
319:         break;
320:       }
321:     }
322:   }
323:   return(0);
324: }

328: /*@
329:   DMLabelSetValue - Set the value a label assigns to a point

331:   Input Parameters:
332: + label - the DMLabel
333: . point - the point
334: - value - The point value

336:   Level: intermediate

338: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue()
339: @*/
340: PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
341: {
342:   PETSC_UNUSED PetscHashIIter iter, ret;
343:   PetscInt                    v;
344:   PetscErrorCode              ierr;

347:   DMLabelMakeInvalid_Private(label);
348:   /* Find, or add, label value */
349:   for (v = 0; v < label->numStrata; ++v) {
350:     if (label->stratumValues[v] == value) break;
351:   }
352:   /* Create new table */
353:   if (v >= label->numStrata) {
354:     PetscInt   *tmpV;
355:     PetscHashI *tmpH;

357:     PetscMalloc1((label->numStrata+1), &tmpV);
358:     PetscMalloc1((label->numStrata+1), &tmpH);
359:     for (v = 0; v < label->numStrata; ++v) {
360:       tmpV[v] = label->stratumValues[v];
361:       tmpH[v] = label->ht[v];
362:     }
363:     tmpV[v] = value;
364:     PetscHashICreate(tmpH[v]);
365:     ++label->numStrata;
366:     PetscFree(label->stratumValues);
367:     PetscFree(label->ht);
368:     label->stratumValues = tmpV;
369:     label->ht            = tmpH;
370:   }
371:   /* Set key */
372:   PetscHashIPut(label->ht[v], point, ret, iter);
373:   return(0);
374: }

378: /*@
379:   DMLabelClearValue - Clear the value a label assigns to a point

381:   Input Parameters:
382: + label - the DMLabel
383: . point - the point
384: - value - The point value

386:   Level: intermediate

388: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue()
389: @*/
390: PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
391: {
392:   PetscInt       v, p;

396:   /* Find label value */
397:   for (v = 0; v < label->numStrata; ++v) {
398:     if (label->stratumValues[v] == value) break;
399:   }
400:   if (v >= label->numStrata) return(0);
401:   if (label->arrayValid) {
402:     /* Check whether point exists */
403:     PetscFindInt(point, label->stratumSizes[v], &label->points[label->stratumOffsets[v]], &p);
404:     if (p >= 0) {
405:       PetscMemmove(&label->points[p+label->stratumOffsets[v]], &label->points[p+label->stratumOffsets[v]+1], (label->stratumSizes[v]-p-1) * sizeof(PetscInt));
406:       --label->stratumSizes[v];
407:       if (label->bt) {
408:         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);
409:         PetscBTClear(label->bt, point - label->pStart);
410:       }
411:     }
412:   } else {
413:     PetscHashIDelKey(label->ht[v], point);
414:   }
415:   return(0);
416: }

420: PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
421: {

426:   DMLabelMakeValid_Private(label);
427:   *numValues = label->numStrata;
428:   return(0);
429: }

433: PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
434: {

439:   DMLabelMakeValid_Private(label);
440:   ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);
441:   return(0);
442: }

446: PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
447: {
448:   PetscInt       v;

453:   DMLabelMakeValid_Private(label);
454:   *size = 0;
455:   for (v = 0; v < label->numStrata; ++v) {
456:     if (label->stratumValues[v] == value) {
457:       *size = label->stratumSizes[v];
458:       break;
459:     }
460:   }
461:   return(0);
462: }

466: PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
467: {
468:   PetscInt       v;

474:   DMLabelMakeValid_Private(label);
475:   for (v = 0; v < label->numStrata; ++v) {
476:     if (label->stratumValues[v] != value) continue;
477:     if (start) *start = label->points[label->stratumOffsets[v]];
478:     if (end)   *end   = label->points[label->stratumOffsets[v]+label->stratumSizes[v]-1]+1;
479:     break;
480:   }
481:   return(0);
482: }

486: PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
487: {
488:   PetscInt       v;

493:   DMLabelMakeValid_Private(label);
494:   *points = NULL;
495:   for (v = 0; v < label->numStrata; ++v) {
496:     if (label->stratumValues[v] == value) {
497:       if (label->arrayValid) {
498:         ISCreateGeneral(PETSC_COMM_SELF, label->stratumSizes[v], &label->points[label->stratumOffsets[v]], PETSC_COPY_VALUES, points);
499:       } else {
500:         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Need to implement this to speedup Stratify");
501:       }
502:       break;
503:     }
504:   }
505:   return(0);
506: }

510: PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
511: {
512:   PetscInt       v;

516:   for (v = 0; v < label->numStrata; ++v) {
517:     if (label->stratumValues[v] == value) break;
518:   }
519:   if (v >= label->numStrata) return(0);
520:   if (label->bt) {
521:     PetscInt i;

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

526:       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);
527:       PetscBTClear(label->bt, point - label->pStart);
528:     }
529:   }
530:   if (label->arrayValid) {
531:     label->stratumSizes[v] = 0;
532:   } else {
533:     PetscHashIClear(label->ht[v]);
534:   }
535:   return(0);
536: }

540: PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
541: {
542:   PetscInt       v;

546:   DMLabelMakeValid_Private(label);
547:   label->pStart = start;
548:   label->pEnd   = end;
549:   if (label->bt) {PetscBTDestroy(&label->bt);}
550:   /* Could squish offsets, but would only make sense if I reallocate the storage */
551:   for (v = 0; v < label->numStrata; ++v) {
552:     const PetscInt offset = label->stratumOffsets[v];
553:     const PetscInt size   = label->stratumSizes[v];
554:     PetscInt       off    = offset, q;

556:     for (q = offset; q < offset+size; ++q) {
557:       const PetscInt point = label->points[q];

559:       if ((point < start) || (point >= end)) continue;
560:       label->points[off++] = point;
561:     }
562:     label->stratumSizes[v] = off-offset;
563:   }
564:   DMLabelCreateIndex(label, start, end);
565:   return(0);
566: }

570: PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
571: {
572:   const PetscInt *perm;
573:   PetscInt        numValues, numPoints, v, q;
574:   PetscErrorCode  ierr;

577:   DMLabelMakeValid_Private(label);
578:   DMLabelDuplicate(label, labelNew);
579:   DMLabelGetNumValues(*labelNew, &numValues);
580:   ISGetLocalSize(permutation, &numPoints);
581:   ISGetIndices(permutation, &perm);
582:   for (v = 0; v < numValues; ++v) {
583:     const PetscInt offset = (*labelNew)->stratumOffsets[v];
584:     const PetscInt size   = (*labelNew)->stratumSizes[v];

586:     for (q = offset; q < offset+size; ++q) {
587:       const PetscInt point = (*labelNew)->points[q];

589:       if ((point < 0) || (point >= numPoints)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %d is not in [0, %d) for the remapping", point, numPoints);
590:       (*labelNew)->points[q] = perm[point];
591:     }
592:     PetscSortInt(size, &(*labelNew)->points[offset]);
593:   }
594:   ISRestoreIndices(permutation, &perm);
595:   if (label->bt) {
596:     PetscBTDestroy(&label->bt);
597:     DMLabelCreateIndex(label, label->pStart, label->pEnd);
598:   }
599:   return(0);
600: }

604: PetscErrorCode DMLabelDistribute(DMLabel label, PetscSection partSection, IS part, ISLocalToGlobalMapping renumbering, DMLabel *labelNew)
605: {
606:   MPI_Comm       comm = PetscObjectComm((PetscObject) partSection);
607:   PetscInt      *stratumSizes = NULL, *points = NULL, s, p;
608:   PetscMPIInt   *sendcnts = NULL, *offsets = NULL, *displs = NULL, proc;
609:   char          *name;
610:   PetscInt       nameSize;
611:   size_t         len = 0;
612:   PetscMPIInt    rank, numProcs;

616:   if (label) {DMLabelMakeValid_Private(label);}
617:   MPI_Comm_rank(comm, &rank);
618:   MPI_Comm_size(comm, &numProcs);
619:   /* Bcast name */
620:   if (!rank) {PetscStrlen(label->name, &len);}
621:   nameSize = len;
622:   MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);
623:   PetscMalloc(nameSize+1, &name);
624:   if (!rank) {PetscMemcpy(name, label->name, nameSize+1);}
625:   MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);
626:   DMLabelCreate(name, labelNew);
627:   PetscFree(name);
628:   /* Bcast numStrata */
629:   if (!rank) (*labelNew)->numStrata = label->numStrata;
630:   MPI_Bcast(&(*labelNew)->numStrata, 1, MPIU_INT, 0, comm);
631:   PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);
632:   PetscMalloc2((*labelNew)->numStrata,&(*labelNew)->stratumSizes,(*labelNew)->numStrata+1,&(*labelNew)->stratumOffsets);
633:   /* Bcast stratumValues */
634:   if (!rank) {PetscMemcpy((*labelNew)->stratumValues, label->stratumValues, label->numStrata * sizeof(PetscInt));}
635:   MPI_Bcast((*labelNew)->stratumValues, (*labelNew)->numStrata, MPIU_INT, 0, comm);
636:   /* Find size on each process and Scatter: we use the fact that both the stratum points and partArray are sorted */
637:   if (!rank) {
638:     const PetscInt *partArray;
639:     PetscInt        proc;

641:     ISGetIndices(part, &partArray);
642:     PetscCalloc1(numProcs*label->numStrata, &stratumSizes);
643:     /* TODO We should switch to using binary search if the label is a lot smaller than partitions */
644:     for (proc = 0; proc < numProcs; ++proc) {
645:       PetscInt dof, off;

647:       PetscSectionGetDof(partSection, proc, &dof);
648:       PetscSectionGetOffset(partSection, proc, &off);
649:       for (s = 0; s < label->numStrata; ++s) {
650:         PetscInt lStart = label->stratumOffsets[s], lEnd = label->stratumOffsets[s]+label->stratumSizes[s];
651:         PetscInt pStart = off,                      pEnd = off+dof;

653:         while (pStart < pEnd && lStart < lEnd) {
654:           if (partArray[pStart] > label->points[lStart]) {
655:             ++lStart;
656:           } else if (label->points[lStart] > partArray[pStart]) {
657:             ++pStart;
658:           } else {
659:             ++stratumSizes[proc*label->numStrata+s];
660:             ++pStart; ++lStart;
661:           }
662:         }
663:       }
664:     }
665:     ISRestoreIndices(part, &partArray);
666:   }
667:   MPI_Scatter(stratumSizes, (*labelNew)->numStrata, MPIU_INT, (*labelNew)->stratumSizes, (*labelNew)->numStrata, MPIU_INT, 0, comm);
668:   /* Calculate stratumOffsets */
669:   (*labelNew)->stratumOffsets[0] = 0;
670:   for (s = 0; s < (*labelNew)->numStrata; ++s) {(*labelNew)->stratumOffsets[s+1] = (*labelNew)->stratumSizes[s] + (*labelNew)->stratumOffsets[s];}
671:   /* Pack points and Scatter */
672:   if (!rank) {
673:     const PetscInt *partArray;

675:     ISGetIndices(part, &partArray);
676:     PetscMalloc3(numProcs,&sendcnts,numProcs,&offsets,numProcs+1,&displs);
677:     displs[0] = 0;
678:     for (p = 0; p < numProcs; ++p) {
679:       sendcnts[p] = 0;
680:       for (s = 0; s < label->numStrata; ++s) sendcnts[p] += stratumSizes[p*label->numStrata+s];
681:       offsets[p]  = displs[p];
682:       displs[p+1] = displs[p] + sendcnts[p];
683:     }
684:     PetscMalloc1(displs[numProcs], &points);
685:     /* TODO We should switch to using binary search if the label is a lot smaller than partitions */
686:     for (proc = 0; proc < numProcs; ++proc) {
687:       PetscInt dof, off;

689:       PetscSectionGetDof(partSection, proc, &dof);
690:       PetscSectionGetOffset(partSection, proc, &off);
691:       for (s = 0; s < label->numStrata; ++s) {
692:         PetscInt lStart = label->stratumOffsets[s], lEnd = label->stratumOffsets[s]+label->stratumSizes[s];
693:         PetscInt pStart = off,                     pEnd = off+dof;

695:         while (pStart < pEnd && lStart < lEnd) {
696:           if (partArray[pStart] > label->points[lStart]) {
697:             ++lStart;
698:           } else if (label->points[lStart] > partArray[pStart]) {
699:             ++pStart;
700:           } else {
701:             points[offsets[proc]++] = label->points[lStart];
702:             ++pStart; ++lStart;
703:           }
704:         }
705:       }
706:     }
707:     ISRestoreIndices(part, &partArray);
708:   }
709:   PetscMalloc1((*labelNew)->stratumOffsets[(*labelNew)->numStrata], &(*labelNew)->points);
710:   MPI_Scatterv(points, sendcnts, displs, MPIU_INT, (*labelNew)->points, (*labelNew)->stratumOffsets[(*labelNew)->numStrata], MPIU_INT, 0, comm);
711:   PetscFree(points);
712:   PetscFree3(sendcnts,offsets,displs);
713:   PetscFree(stratumSizes);
714:   /* Renumber points */
715:   ISGlobalToLocalMappingApply(renumbering, IS_GTOLM_MASK, (*labelNew)->stratumOffsets[(*labelNew)->numStrata], (*labelNew)->points, NULL, (*labelNew)->points);
716:   /* Sort points */
717:   for (s = 0; s < (*labelNew)->numStrata; ++s) {PetscSortInt((*labelNew)->stratumSizes[s], &(*labelNew)->points[(*labelNew)->stratumOffsets[s]]);}
718:   return(0);
719: }


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

727:   Not Collective

729:   Input Parameters:
730: + dm   - The DMPlex object
731: - name - The label name

733:   Level: intermediate

735: .keywords: mesh
736: .seealso: DMLabelCreate(), DMPlexHasLabel(), DMPlexGetLabelValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
737: @*/
738: PetscErrorCode DMPlexCreateLabel(DM dm, const char name[])
739: {
740:   DM_Plex        *mesh = (DM_Plex*) dm->data;
741:   DMLabel        next  = mesh->labels;
742:   PetscBool      flg   = PETSC_FALSE;

748:   while (next) {
749:     PetscStrcmp(name, next->name, &flg);
750:     if (flg) break;
751:     next = next->next;
752:   }
753:   if (!flg) {
754:     DMLabel tmpLabel = mesh->labels;

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

758:     mesh->labels->next = tmpLabel;
759:   }
760:   return(0);
761: }

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

768:   Not Collective

770:   Input Parameters:
771: + dm   - The DMPlex object
772: . name - The label name
773: - point - The mesh point

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

778:   Level: beginner

780: .keywords: mesh
781: .seealso: DMLabelGetValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
782: @*/
783: PetscErrorCode DMPlexGetLabelValue(DM dm, const char name[], PetscInt point, PetscInt *value)
784: {
785:   DMLabel        label;

791:   DMPlexGetLabel(dm, name, &label);
792:   if (!label) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No label named %s was found", name);
793:   DMLabelGetValue(label, point, value);
794:   return(0);
795: }

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

802:   Not Collective

804:   Input Parameters:
805: + dm   - The DMPlex object
806: . name - The label name
807: . point - The mesh point
808: - value - The label value for this point

810:   Output Parameter:

812:   Level: beginner

814: .keywords: mesh
815: .seealso: DMLabelSetValue(), DMPlexGetStratumIS(), DMPlexClearLabelValue()
816: @*/
817: PetscErrorCode DMPlexSetLabelValue(DM dm, const char name[], PetscInt point, PetscInt value)
818: {
819:   DMLabel        label;

825:   DMPlexGetLabel(dm, name, &label);
826:   if (!label) {
827:     DMPlexCreateLabel(dm, name);
828:     DMPlexGetLabel(dm, name, &label);
829:   }
830:   DMLabelSetValue(label, point, value);
831:   return(0);
832: }

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

839:   Not Collective

841:   Input Parameters:
842: + dm   - The DMPlex object
843: . name - The label name
844: . point - The mesh point
845: - value - The label value for this point

847:   Output Parameter:

849:   Level: beginner

851: .keywords: mesh
852: .seealso: DMLabelClearValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
853: @*/
854: PetscErrorCode DMPlexClearLabelValue(DM dm, const char name[], PetscInt point, PetscInt value)
855: {
856:   DMLabel        label;

862:   DMPlexGetLabel(dm, name, &label);
863:   if (!label) return(0);
864:   DMLabelClearValue(label, point, value);
865:   return(0);
866: }

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

873:   Not Collective

875:   Input Parameters:
876: + dm   - The DMPlex object
877: - name - The label name

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

882:   Level: beginner

884: .keywords: mesh
885: .seealso: DMLabeGetNumValues(), DMPlexSetLabelValue()
886: @*/
887: PetscErrorCode DMPlexGetLabelSize(DM dm, const char name[], PetscInt *size)
888: {
889:   DMLabel        label;

896:   DMPlexGetLabel(dm, name, &label);
897:   *size = 0;
898:   if (!label) return(0);
899:   DMLabelGetNumValues(label, size);
900:   return(0);
901: }

905: /*@C
906:   DMPlexGetLabelIdIS - Get the integer ids in a label

908:   Not Collective

910:   Input Parameters:
911: + mesh - The DMPlex object
912: - name - The label name

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

917:   Level: beginner

919: .keywords: mesh
920: .seealso: DMLabelGetValueIS(), DMPlexGetLabelSize()
921: @*/
922: PetscErrorCode DMPlexGetLabelIdIS(DM dm, const char name[], IS *ids)
923: {
924:   DMLabel        label;

931:   DMPlexGetLabel(dm, name, &label);
932:   *ids = NULL;
933:   if (!label) return(0);
934:   DMLabelGetValueIS(label, ids);
935:   return(0);
936: }

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

943:   Not Collective

945:   Input Parameters:
946: + dm - The DMPlex object
947: . name - The label name
948: - value - The stratum value

950:   Output Parameter:
951: . size - The stratum size

953:   Level: beginner

955: .keywords: mesh
956: .seealso: DMLabelGetStratumSize(), DMPlexGetLabelSize(), DMPlexGetLabelIds()
957: @*/
958: PetscErrorCode DMPlexGetStratumSize(DM dm, const char name[], PetscInt value, PetscInt *size)
959: {
960:   DMLabel        label;

967:   DMPlexGetLabel(dm, name, &label);
968:   *size = 0;
969:   if (!label) return(0);
970:   DMLabelGetStratumSize(label, value, size);
971:   return(0);
972: }

976: /*@C
977:   DMPlexGetStratumIS - Get the points in a label stratum

979:   Not Collective

981:   Input Parameters:
982: + dm - The DMPlex object
983: . name - The label name
984: - value - The stratum value

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

989:   Level: beginner

991: .keywords: mesh
992: .seealso: DMLabelGetStratumIS(), DMPlexGetStratumSize()
993: @*/
994: PetscErrorCode DMPlexGetStratumIS(DM dm, const char name[], PetscInt value, IS *points)
995: {
996:   DMLabel        label;

1003:   DMPlexGetLabel(dm, name, &label);
1004:   *points = NULL;
1005:   if (!label) return(0);
1006:   DMLabelGetStratumIS(label, value, points);
1007:   return(0);
1008: }

1012: /*@C
1013:   DMPlexClearLabelStratum - Remove all points from a stratum from a Sieve Label

1015:   Not Collective

1017:   Input Parameters:
1018: + dm   - The DMPlex object
1019: . name - The label name
1020: - value - The label value for this point

1022:   Output Parameter:

1024:   Level: beginner

1026: .keywords: mesh
1027: .seealso: DMLabelClearStratum(), DMPlexSetLabelValue(), DMPlexGetStratumIS(), DMPlexClearLabelValue()
1028: @*/
1029: PetscErrorCode DMPlexClearLabelStratum(DM dm, const char name[], PetscInt value)
1030: {
1031:   DMLabel        label;

1037:   DMPlexGetLabel(dm, name, &label);
1038:   if (!label) return(0);
1039:   DMLabelClearStratum(label, value);
1040:   return(0);
1041: }

1045: /*@
1046:   DMPlexGetNumLabels - Return the number of labels defined by the mesh

1048:   Not Collective

1050:   Input Parameter:
1051: . dm   - The DMPlex object

1053:   Output Parameter:
1054: . numLabels - the number of Labels

1056:   Level: intermediate

1058: .keywords: mesh
1059: .seealso: DMPlexGetLabelValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
1060: @*/
1061: PetscErrorCode DMPlexGetNumLabels(DM dm, PetscInt *numLabels)
1062: {
1063:   DM_Plex  *mesh = (DM_Plex*) dm->data;
1064:   DMLabel  next  = mesh->labels;
1065:   PetscInt n     = 0;

1070:   while (next) {
1071:     ++n;
1072:     next = next->next;
1073:   }
1074:   *numLabels = n;
1075:   return(0);
1076: }

1080: /*@C
1081:   DMPlexGetLabelName - Return the name of nth label

1083:   Not Collective

1085:   Input Parameters:
1086: + dm - The DMPlex object
1087: - n  - the label number

1089:   Output Parameter:
1090: . name - the label name

1092:   Level: intermediate

1094: .keywords: mesh
1095: .seealso: DMPlexGetLabelValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
1096: @*/
1097: PetscErrorCode DMPlexGetLabelName(DM dm, PetscInt n, const char **name)
1098: {
1099:   DM_Plex  *mesh = (DM_Plex*) dm->data;
1100:   DMLabel  next  = mesh->labels;
1101:   PetscInt l     = 0;

1106:   while (next) {
1107:     if (l == n) {
1108:       *name = next->name;
1109:       return(0);
1110:     }
1111:     ++l;
1112:     next = next->next;
1113:   }
1114:   SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %d does not exist in this DM", n);
1115: }

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

1122:   Not Collective

1124:   Input Parameters:
1125: + dm   - The DMPlex object
1126: - name - The label name

1128:   Output Parameter:
1129: . hasLabel - PETSC_TRUE if the label is present

1131:   Level: intermediate

1133: .keywords: mesh
1134: .seealso: DMPlexCreateLabel(), DMPlexGetLabelValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
1135: @*/
1136: PetscErrorCode DMPlexHasLabel(DM dm, const char name[], PetscBool *hasLabel)
1137: {
1138:   DM_Plex        *mesh = (DM_Plex*) dm->data;
1139:   DMLabel        next  = mesh->labels;

1146:   *hasLabel = PETSC_FALSE;
1147:   while (next) {
1148:     PetscStrcmp(name, next->name, hasLabel);
1149:     if (*hasLabel) break;
1150:     next = next->next;
1151:   }
1152:   return(0);
1153: }

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

1160:   Not Collective

1162:   Input Parameters:
1163: + dm   - The DMPlex object
1164: - name - The label name

1166:   Output Parameter:
1167: . label - The DMLabel, or NULL if the label is absent

1169:   Level: intermediate

1171: .keywords: mesh
1172: .seealso: DMPlexCreateLabel(), DMPlexHasLabel(), DMPlexGetLabelValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
1173: @*/
1174: PetscErrorCode DMPlexGetLabel(DM dm, const char name[], DMLabel *label)
1175: {
1176:   DM_Plex        *mesh = (DM_Plex*) dm->data;
1177:   DMLabel        next  = mesh->labels;
1178:   PetscBool      hasLabel;

1185:   *label = NULL;
1186:   while (next) {
1187:     PetscStrcmp(name, next->name, &hasLabel);
1188:     if (hasLabel) {
1189:       *label = next;
1190:       break;
1191:     }
1192:     next = next->next;
1193:   }
1194:   return(0);
1195: }

1199: /*@C
1200:   DMPlexAddLabel - Add the label to this mesh

1202:   Not Collective

1204:   Input Parameters:
1205: + dm   - The DMPlex object
1206: - label - The DMLabel

1208:   Level: developer

1210: .keywords: mesh
1211: .seealso: DMPlexCreateLabel(), DMPlexHasLabel(), DMPlexGetLabelValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
1212: @*/
1213: PetscErrorCode DMPlexAddLabel(DM dm, DMLabel label)
1214: {
1215:   DM_Plex        *mesh = (DM_Plex*) dm->data;
1216:   PetscBool      hasLabel;

1221:   DMPlexHasLabel(dm, label->name, &hasLabel);
1222:   if (hasLabel) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %s already exists in this DM", label->name);
1223:   label->next  = mesh->labels;
1224:   mesh->labels = label;
1225:   return(0);
1226: }

1230: /*@C
1231:   DMPlexRemoveLabel - Remove the label from this mesh

1233:   Not Collective

1235:   Input Parameters:
1236: + dm   - The DMPlex object
1237: - name - The label name

1239:   Output Parameter:
1240: . label - The DMLabel, or NULL if the label is absent

1242:   Level: developer

1244: .keywords: mesh
1245: .seealso: DMPlexCreateLabel(), DMPlexHasLabel(), DMPlexGetLabelValue(), DMPlexSetLabelValue(), DMPlexGetStratumIS()
1246: @*/
1247: PetscErrorCode DMPlexRemoveLabel(DM dm, const char name[], DMLabel *label)
1248: {
1249:   DM_Plex        *mesh = (DM_Plex*) dm->data;
1250:   DMLabel        next  = mesh->labels;
1251:   DMLabel        last  = NULL;
1252:   PetscBool      hasLabel;

1257:   DMPlexHasLabel(dm, name, &hasLabel);
1258:   *label = NULL;
1259:   if (!hasLabel) return(0);
1260:   while (next) {
1261:     PetscStrcmp(name, next->name, &hasLabel);
1262:     if (hasLabel) {
1263:       if (last) last->next   = next->next;
1264:       else      mesh->labels = next->next;
1265:       next->next = NULL;
1266:       *label     = next;
1267:       break;
1268:     }
1269:     last = next;
1270:     next = next->next;
1271:   }
1272:   return(0);
1273: }