Actual source code: dmlabel.c
petsc-3.11.0 2019-03-29
1: #include <petscdm.h>
2: #include <petsc/private/dmlabelimpl.h>
3: #include <petsc/private/isimpl.h>
4: #include <petscsf.h>
6: /*@C
7: DMLabelCreate - Create a DMLabel object, which is a multimap
9: Input parameters:
10: + comm - The communicator, usually PETSC_COMM_SELF
11: - name - The label name
13: Output parameter:
14: . label - The DMLabel
16: Level: beginner
18: .seealso: DMLabelDestroy()
19: @*/
20: PetscErrorCode DMLabelCreate(MPI_Comm comm, const char name[], DMLabel *label)
21: {
26: DMInitializePackage();
28: PetscHeaderCreate(*label,DMLABEL_CLASSID,"DMLabel","DMLabel","DM",comm,DMLabelDestroy,DMLabelView);
30: (*label)->numStrata = 0;
31: (*label)->defaultValue = -1;
32: (*label)->stratumValues = NULL;
33: (*label)->validIS = NULL;
34: (*label)->stratumSizes = NULL;
35: (*label)->points = NULL;
36: (*label)->ht = NULL;
37: (*label)->pStart = -1;
38: (*label)->pEnd = -1;
39: (*label)->bt = NULL;
40: PetscHMapICreate(&(*label)->hmap);
41: PetscObjectSetName((PetscObject) *label, name);
42: return(0);
43: }
45: /*
46: DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format
48: Input parameter:
49: + label - The DMLabel
50: - v - The stratum value
52: Output parameter:
53: . label - The DMLabel with stratum in sorted list format
55: Level: developer
57: .seealso: DMLabelCreate()
58: */
59: static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v)
60: {
61: PetscInt off = 0, *pointArray, p;
64: if (PetscLikely(v >= 0 && v < label->numStrata) && label->validIS[v]) return 0;
66: if (v < 0 || v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeValid_Private\n", v);
67: PetscHSetIGetSize(label->ht[v], &label->stratumSizes[v]);
68: PetscMalloc1(label->stratumSizes[v], &pointArray);
69: PetscHSetIGetElems(label->ht[v], &off, pointArray);
70: PetscHSetIClear(label->ht[v]);
71: PetscSortInt(label->stratumSizes[v], pointArray);
72: if (label->bt) {
73: for (p = 0; p < label->stratumSizes[v]; ++p) {
74: const PetscInt point = pointArray[p];
75: 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);
76: PetscBTSet(label->bt, point - label->pStart);
77: }
78: }
79: ISCreateGeneral(PETSC_COMM_SELF,label->stratumSizes[v],pointArray,PETSC_OWN_POINTER,&(label->points[v]));
80: PetscObjectSetName((PetscObject) (label->points[v]), "indices");
81: label->validIS[v] = PETSC_TRUE;
82: PetscObjectStateIncrease((PetscObject) label);
83: return(0);
84: }
86: /*
87: DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format
89: Input parameter:
90: . label - The DMLabel
92: Output parameter:
93: . label - The DMLabel with all strata in sorted list format
95: Level: developer
97: .seealso: DMLabelCreate()
98: */
99: static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label)
100: {
101: PetscInt v;
105: for (v = 0; v < label->numStrata; v++){
106: DMLabelMakeValid_Private(label, v);
107: }
108: return(0);
109: }
111: /*
112: DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format
114: Input parameter:
115: + label - The DMLabel
116: - v - The stratum value
118: Output parameter:
119: . label - The DMLabel with stratum in hash format
121: Level: developer
123: .seealso: DMLabelCreate()
124: */
125: static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v)
126: {
127: PetscInt p;
128: const PetscInt *points;
131: if (PetscLikely(v >= 0 && v < label->numStrata) && !label->validIS[v]) return 0;
133: if (v < 0 || v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeInvalid_Private\n", v);
134: if (label->points[v]) {
135: ISGetIndices(label->points[v], &points);
136: for (p = 0; p < label->stratumSizes[v]; ++p) {
137: PetscHSetIAdd(label->ht[v], points[p]);
138: }
139: ISRestoreIndices(label->points[v],&points);
140: ISDestroy(&(label->points[v]));
141: }
142: label->validIS[v] = PETSC_FALSE;
143: return(0);
144: }
146: #if !defined(DMLABEL_LOOKUP_THRESHOLD)
147: #define DMLABEL_LOOKUP_THRESHOLD 16
148: #endif
150: PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupStratum(DMLabel label, PetscInt value, PetscInt *index)
151: {
152: PetscInt v;
156: *index = -1;
157: if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD) {
158: for (v = 0; v < label->numStrata; ++v)
159: if (label->stratumValues[v] == value) {*index = v; break;}
160: } else {
161: PetscHMapIGet(label->hmap, value, index);
162: #if defined(PETSC_USE_DEBUG)
163: v = *index;
164: if (v >= 0 && v >= label->numStrata) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent lookup using strata hash map");
165: if (v >= 0 && label->stratumValues[v] != value) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent lookup using strata hash map");
166: #endif
167: }
168: return(0);
169: }
171: PETSC_STATIC_INLINE PetscErrorCode DMLabelNewStratum(DMLabel label, PetscInt value, PetscInt *index)
172: {
173: PetscInt v;
174: PetscInt *tmpV;
175: PetscInt *tmpS;
176: PetscHSetI *tmpH, ht;
177: IS *tmpP, is;
178: PetscBool *tmpB;
179: PetscHMapI hmap = label->hmap;
183: v = label->numStrata;
184: tmpV = label->stratumValues;
185: tmpS = label->stratumSizes;
186: tmpH = label->ht;
187: tmpP = label->points;
188: tmpB = label->validIS;
189: { /* TODO: PetscRealloc() is broken, use malloc+memcpy+free */
190: PetscInt *oldV = tmpV;
191: PetscInt *oldS = tmpS;
192: PetscHSetI *oldH = tmpH;
193: IS *oldP = tmpP;
194: PetscBool *oldB = tmpB;
195: PetscMalloc((v+1)*sizeof(*tmpV), &tmpV);
196: PetscMalloc((v+1)*sizeof(*tmpS), &tmpS);
197: PetscMalloc((v+1)*sizeof(*tmpH), &tmpH);
198: PetscMalloc((v+1)*sizeof(*tmpP), &tmpP);
199: PetscMalloc((v+1)*sizeof(*tmpB), &tmpB);
200: PetscMemcpy(tmpV, oldV, v*sizeof(*tmpV));
201: PetscMemcpy(tmpS, oldS, v*sizeof(*tmpS));
202: PetscMemcpy(tmpH, oldH, v*sizeof(*tmpH));
203: PetscMemcpy(tmpP, oldP, v*sizeof(*tmpP));
204: PetscMemcpy(tmpB, oldB, v*sizeof(*tmpB));
205: PetscFree(oldV);
206: PetscFree(oldS);
207: PetscFree(oldH);
208: PetscFree(oldP);
209: PetscFree(oldB);
210: }
211: label->numStrata = v+1;
212: label->stratumValues = tmpV;
213: label->stratumSizes = tmpS;
214: label->ht = tmpH;
215: label->points = tmpP;
216: label->validIS = tmpB;
217: PetscHSetICreate(&ht);
218: ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);
219: PetscHMapISet(hmap, value, v);
220: tmpV[v] = value;
221: tmpS[v] = 0;
222: tmpH[v] = ht;
223: tmpP[v] = is;
224: tmpB[v] = PETSC_TRUE;
225: *index = v;
226: return(0);
227: }
229: PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupAddStratum(DMLabel label, PetscInt value, PetscInt *index)
230: {
233: DMLabelLookupStratum(label, value, index);
234: if (*index < 0) {DMLabelNewStratum(label, value, index);}
235: return(0);
236: }
238: /*@
239: DMLabelAddStratum - Adds a new stratum value in a DMLabel
241: Input Parameter:
242: + label - The DMLabel
243: - value - The stratum value
245: Level: beginner
247: .seealso: DMLabelCreate(), DMLabelDestroy()
248: @*/
249: PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value)
250: {
251: PetscInt v;
256: DMLabelLookupAddStratum(label, value, &v);
257: return(0);
258: }
260: /*@
261: DMLabelAddStrata - Adds new stratum values in a DMLabel
263: Input Parameter:
264: + label - The DMLabel
265: . numStrata - The number of stratum values
266: - stratumValues - The stratum values
268: Level: beginner
270: .seealso: DMLabelCreate(), DMLabelDestroy()
271: @*/
272: PetscErrorCode DMLabelAddStrata(DMLabel label, PetscInt numStrata, const PetscInt stratumValues[])
273: {
274: PetscInt *values, v;
280: PetscMalloc1(numStrata, &values);
281: PetscMemcpy(values, stratumValues, numStrata*sizeof(PetscInt));
282: PetscSortRemoveDupsInt(&numStrata, values);
283: if (!label->numStrata) { /* Fast preallocation */
284: PetscInt *tmpV;
285: PetscInt *tmpS;
286: PetscHSetI *tmpH, ht;
287: IS *tmpP, is;
288: PetscBool *tmpB;
289: PetscHMapI hmap = label->hmap;
291: PetscMalloc1(numStrata, &tmpV);
292: PetscMalloc1(numStrata, &tmpS);
293: PetscMalloc1(numStrata, &tmpH);
294: PetscMalloc1(numStrata, &tmpP);
295: PetscMalloc1(numStrata, &tmpB);
296: label->numStrata = numStrata;
297: label->stratumValues = tmpV;
298: label->stratumSizes = tmpS;
299: label->ht = tmpH;
300: label->points = tmpP;
301: label->validIS = tmpB;
302: for (v = 0; v < numStrata; ++v) {
303: PetscHSetICreate(&ht);
304: ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);
305: PetscHMapISet(hmap, values[v], v);
306: tmpV[v] = values[v];
307: tmpS[v] = 0;
308: tmpH[v] = ht;
309: tmpP[v] = is;
310: tmpB[v] = PETSC_TRUE;
311: }
312: } else {
313: for (v = 0; v < numStrata; ++v) {
314: DMLabelAddStratum(label, values[v]);
315: }
316: }
317: PetscFree(values);
318: return(0);
319: }
321: /*@
322: DMLabelAddStrataIS - Adds new stratum values in a DMLabel
324: Input Parameter:
325: + label - The DMLabel
326: - valueIS - Index set with stratum values
328: Level: beginner
330: .seealso: DMLabelCreate(), DMLabelDestroy()
331: @*/
332: PetscErrorCode DMLabelAddStrataIS(DMLabel label, IS valueIS)
333: {
334: PetscInt numStrata;
335: const PetscInt *stratumValues;
341: ISGetLocalSize(valueIS, &numStrata);
342: ISGetIndices(valueIS, &stratumValues);
343: DMLabelAddStrata(label, numStrata, stratumValues);
344: return(0);
345: }
347: static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer)
348: {
349: PetscInt v;
350: PetscMPIInt rank;
354: MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);
355: PetscViewerASCIIPushSynchronized(viewer);
356: if (label) {
357: const char *name;
359: PetscObjectGetName((PetscObject) label, &name);
360: PetscViewerASCIIPrintf(viewer, "Label '%s':\n", name);
361: if (label->bt) {PetscViewerASCIIPrintf(viewer, " Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);}
362: for (v = 0; v < label->numStrata; ++v) {
363: const PetscInt value = label->stratumValues[v];
364: const PetscInt *points;
365: PetscInt p;
367: ISGetIndices(label->points[v], &points);
368: for (p = 0; p < label->stratumSizes[v]; ++p) {
369: PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);
370: }
371: ISRestoreIndices(label->points[v],&points);
372: }
373: }
374: PetscViewerFlush(viewer);
375: PetscViewerASCIIPopSynchronized(viewer);
376: return(0);
377: }
379: /*@C
380: DMLabelView - View the label
382: Input Parameters:
383: + label - The DMLabel
384: - viewer - The PetscViewer
386: Level: intermediate
388: .seealso: DMLabelCreate(), DMLabelDestroy()
389: @*/
390: PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
391: {
392: PetscBool iascii;
397: if (!viewer) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)label), &viewer);}
399: if (label) {DMLabelMakeAllValid_Private(label);}
400: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
401: if (iascii) {
402: DMLabelView_Ascii(label, viewer);
403: }
404: return(0);
405: }
407: /*@
408: DMLabelReset - Destroys internal data structures in a DMLabel
410: Input Parameter:
411: . label - The DMLabel
413: Level: beginner
415: .seealso: DMLabelDestroy(), DMLabelCreate()
416: @*/
417: PetscErrorCode DMLabelReset(DMLabel label)
418: {
419: PetscInt v;
424: for (v = 0; v < label->numStrata; ++v) {
425: PetscHSetIDestroy(&label->ht[v]);
426: ISDestroy(&label->points[v]);
427: }
428: label->numStrata = 0;
429: PetscFree(label->stratumValues);
430: PetscFree(label->stratumSizes);
431: PetscFree(label->ht);
432: PetscFree(label->points);
433: PetscFree(label->validIS);
434: PetscHMapIReset(label->hmap);
435: label->pStart = -1;
436: label->pEnd = -1;
437: PetscBTDestroy(&label->bt);
438: return(0);
439: }
441: /*@
442: DMLabelDestroy - Destroys a DMLabel
444: Input Parameter:
445: . label - The DMLabel
447: Level: beginner
449: .seealso: DMLabelReset(), DMLabelCreate()
450: @*/
451: PetscErrorCode DMLabelDestroy(DMLabel *label)
452: {
456: if (!*label) return(0);
458: if (--((PetscObject)(*label))->refct > 0) {*label = NULL; return(0);}
459: DMLabelReset(*label);
460: PetscHMapIDestroy(&(*label)->hmap);
461: PetscHeaderDestroy(label);
462: return(0);
463: }
465: /*@
466: DMLabelDuplicate - Duplicates a DMLabel
468: Input Parameter:
469: . label - The DMLabel
471: Output Parameter:
472: . labelnew - location to put new vector
474: Level: intermediate
476: .seealso: DMLabelCreate(), DMLabelDestroy()
477: @*/
478: PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
479: {
480: const char *name;
481: PetscInt v;
486: DMLabelMakeAllValid_Private(label);
487: PetscObjectGetName((PetscObject) label, &name);
488: DMLabelCreate(PetscObjectComm((PetscObject) label), name, labelnew);
490: (*labelnew)->numStrata = label->numStrata;
491: (*labelnew)->defaultValue = label->defaultValue;
492: PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);
493: PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);
494: PetscMalloc1(label->numStrata, &(*labelnew)->ht);
495: PetscMalloc1(label->numStrata, &(*labelnew)->points);
496: PetscMalloc1(label->numStrata, &(*labelnew)->validIS);
497: for (v = 0; v < label->numStrata; ++v) {
498: PetscHSetICreate(&(*labelnew)->ht[v]);
499: (*labelnew)->stratumValues[v] = label->stratumValues[v];
500: (*labelnew)->stratumSizes[v] = label->stratumSizes[v];
501: PetscObjectReference((PetscObject) (label->points[v]));
502: (*labelnew)->points[v] = label->points[v];
503: (*labelnew)->validIS[v] = PETSC_TRUE;
504: }
505: PetscHMapIDuplicate(label->hmap,&(*labelnew)->hmap);
506: (*labelnew)->pStart = -1;
507: (*labelnew)->pEnd = -1;
508: (*labelnew)->bt = NULL;
509: return(0);
510: }
512: /*@
513: DMLabelComputeIndex - Create an index structure for membership determination, automatically determining the bounds
515: Input Parameter:
516: . label - The DMLabel
518: Level: intermediate
520: .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
521: @*/
522: PetscErrorCode DMLabelComputeIndex(DMLabel label)
523: {
524: PetscInt pStart = PETSC_MAX_INT, pEnd = -1, v;
529: DMLabelMakeAllValid_Private(label);
530: for (v = 0; v < label->numStrata; ++v) {
531: const PetscInt *points;
532: PetscInt i;
534: ISGetIndices(label->points[v], &points);
535: for (i = 0; i < label->stratumSizes[v]; ++i) {
536: const PetscInt point = points[i];
538: pStart = PetscMin(point, pStart);
539: pEnd = PetscMax(point+1, pEnd);
540: }
541: ISRestoreIndices(label->points[v], &points);
542: }
543: label->pStart = pStart == PETSC_MAX_INT ? -1 : pStart;
544: label->pEnd = pEnd;
545: DMLabelCreateIndex(label, label->pStart, label->pEnd);
546: return(0);
547: }
549: /*@
550: DMLabelCreateIndex - Create an index structure for membership determination
552: Input Parameters:
553: + label - The DMLabel
554: . pStart - The smallest point
555: - pEnd - The largest point + 1
557: Level: intermediate
559: .seealso: DMLabelHasPoint(), DMLabelComputeIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
560: @*/
561: PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
562: {
563: PetscInt v;
568: DMLabelDestroyIndex(label);
569: DMLabelMakeAllValid_Private(label);
570: label->pStart = pStart;
571: label->pEnd = pEnd;
572: /* This can be hooked into SetValue(), ClearValue(), etc. for updating */
573: PetscBTCreate(pEnd - pStart, &label->bt);
574: for (v = 0; v < label->numStrata; ++v) {
575: const PetscInt *points;
576: PetscInt i;
578: ISGetIndices(label->points[v], &points);
579: for (i = 0; i < label->stratumSizes[v]; ++i) {
580: const PetscInt point = points[i];
582: if ((point < pStart) || (point >= pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, pStart, pEnd);
583: PetscBTSet(label->bt, point - pStart);
584: }
585: ISRestoreIndices(label->points[v], &points);
586: }
587: return(0);
588: }
590: /*@
591: DMLabelDestroyIndex - Destroy the index structure
593: Input Parameter:
594: . label - the DMLabel
596: Level: intermediate
598: .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
599: @*/
600: PetscErrorCode DMLabelDestroyIndex(DMLabel label)
601: {
606: label->pStart = -1;
607: label->pEnd = -1;
608: PetscBTDestroy(&label->bt);
609: return(0);
610: }
612: /*@
613: DMLabelGetBounds - Return the smallest and largest point in the label
615: Input Parameter:
616: . label - the DMLabel
618: Output Parameters:
619: + pStart - The smallest point
620: - pEnd - The largest point + 1
622: Note: This will compute an index for the label if one does not exist.
624: Level: intermediate
626: .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
627: @*/
628: PetscErrorCode DMLabelGetBounds(DMLabel label, PetscInt *pStart, PetscInt *pEnd)
629: {
634: if ((label->pStart == -1) && (label->pEnd == -1)) {DMLabelComputeIndex(label);}
635: if (pStart) {
637: *pStart = label->pStart;
638: }
639: if (pEnd) {
641: *pEnd = label->pEnd;
642: }
643: return(0);
644: }
646: /*@
647: DMLabelHasValue - Determine whether a label assigns the value to any point
649: Input Parameters:
650: + label - the DMLabel
651: - value - the value
653: Output Parameter:
654: . contains - Flag indicating whether the label maps this value to any point
656: Level: developer
658: .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue()
659: @*/
660: PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains)
661: {
662: PetscInt v;
668: DMLabelLookupStratum(label, value, &v);
669: *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE;
670: return(0);
671: }
673: /*@
674: DMLabelHasPoint - Determine whether a label assigns a value to a point
676: Input Parameters:
677: + label - the DMLabel
678: - point - the point
680: Output Parameter:
681: . contains - Flag indicating whether the label maps this point to a value
683: Note: The user must call DMLabelCreateIndex() before this function.
685: Level: developer
687: .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
688: @*/
689: PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
690: {
696: DMLabelMakeAllValid_Private(label);
697: #if defined(PETSC_USE_DEBUG)
698: if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
699: 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);
700: #endif
701: *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
702: return(0);
703: }
705: /*@
706: DMLabelStratumHasPoint - Return true if the stratum contains a point
708: Input Parameters:
709: + label - the DMLabel
710: . value - the stratum value
711: - point - the point
713: Output Parameter:
714: . contains - true if the stratum contains the point
716: Level: intermediate
718: .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue()
719: @*/
720: PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains)
721: {
722: PetscInt v;
728: *contains = PETSC_FALSE;
729: DMLabelLookupStratum(label, value, &v);
730: if (v < 0) return(0);
732: if (label->validIS[v]) {
733: PetscInt i;
735: ISLocate(label->points[v], point, &i);
736: if (i >= 0) *contains = PETSC_TRUE;
737: } else {
738: PetscBool has;
740: PetscHSetIHas(label->ht[v], point, &has);
741: if (has) *contains = PETSC_TRUE;
742: }
743: return(0);
744: }
746: /*@
747: DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
748: When a label is created, it is initialized to -1.
750: Input parameter:
751: . label - a DMLabel object
753: Output parameter:
754: . defaultValue - the default value
756: Level: beginner
758: .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
759: @*/
760: PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue)
761: {
764: *defaultValue = label->defaultValue;
765: return(0);
766: }
768: /*@
769: DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
770: When a label is created, it is initialized to -1.
772: Input parameter:
773: . label - a DMLabel object
775: Output parameter:
776: . defaultValue - the default value
778: Level: beginner
780: .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
781: @*/
782: PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue)
783: {
786: label->defaultValue = defaultValue;
787: return(0);
788: }
790: /*@
791: DMLabelGetValue - Return the value a label assigns to a point, or the label's default value (which is initially -1, and can be changed with DMLabelSetDefaultValue())
793: Input Parameters:
794: + label - the DMLabel
795: - point - the point
797: Output Parameter:
798: . value - The point value, or the default value (-1 by default)
800: Level: intermediate
802: .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
803: @*/
804: PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
805: {
806: PetscInt v;
812: *value = label->defaultValue;
813: for (v = 0; v < label->numStrata; ++v) {
814: if (label->validIS[v]) {
815: PetscInt i;
817: ISLocate(label->points[v], point, &i);
818: if (i >= 0) {
819: *value = label->stratumValues[v];
820: break;
821: }
822: } else {
823: PetscBool has;
825: PetscHSetIHas(label->ht[v], point, &has);
826: if (has) {
827: *value = label->stratumValues[v];
828: break;
829: }
830: }
831: }
832: return(0);
833: }
835: /*@
836: DMLabelSetValue - Set the value a label assigns to a point. If the value is the same as the label's default value (which is initially -1, and can be changed with DMLabelSetDefaultValue() to something different), then this function will do nothing.
838: Input Parameters:
839: + label - the DMLabel
840: . point - the point
841: - value - The point value
843: Level: intermediate
845: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
846: @*/
847: PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
848: {
849: PetscInt v;
854: /* Find label value, add new entry if needed */
855: if (value == label->defaultValue) return(0);
856: DMLabelLookupAddStratum(label, value, &v);
857: /* Set key */
858: DMLabelMakeInvalid_Private(label, v);
859: PetscHSetIAdd(label->ht[v], point);
860: return(0);
861: }
863: /*@
864: DMLabelClearValue - Clear the value a label assigns to a point
866: Input Parameters:
867: + label - the DMLabel
868: . point - the point
869: - value - The point value
871: Level: intermediate
873: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue()
874: @*/
875: PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
876: {
877: PetscInt v;
882: /* Find label value */
883: DMLabelLookupStratum(label, value, &v);
884: if (v < 0) return(0);
886: if (label->bt) {
887: 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);
888: PetscBTClear(label->bt, point - label->pStart);
889: }
891: /* Delete key */
892: DMLabelMakeInvalid_Private(label, v);
893: PetscHSetIDel(label->ht[v], point);
894: return(0);
895: }
897: /*@
898: DMLabelInsertIS - Set all points in the IS to a value
900: Input Parameters:
901: + label - the DMLabel
902: . is - the point IS
903: - value - The point value
905: Level: intermediate
907: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
908: @*/
909: PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value)
910: {
911: PetscInt v, n, p;
912: const PetscInt *points;
913: PetscErrorCode ierr;
918: /* Find label value, add new entry if needed */
919: if (value == label->defaultValue) return(0);
920: DMLabelLookupAddStratum(label, value, &v);
921: /* Set keys */
922: DMLabelMakeInvalid_Private(label, v);
923: ISGetLocalSize(is, &n);
924: ISGetIndices(is, &points);
925: for (p = 0; p < n; ++p) {PetscHSetIAdd(label->ht[v], points[p]);}
926: ISRestoreIndices(is, &points);
927: return(0);
928: }
930: /*@
931: DMLabelGetNumValues - Get the number of values that the DMLabel takes
933: Input Parameter:
934: . label - the DMLabel
936: Output Paramater:
937: . numValues - the number of values
939: Level: intermediate
941: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
942: @*/
943: PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
944: {
948: *numValues = label->numStrata;
949: return(0);
950: }
952: /*@
953: DMLabelGetValueIS - Get an IS of all values that the DMlabel takes
955: Input Parameter:
956: . label - the DMLabel
958: Output Paramater:
959: . is - the value IS
961: Level: intermediate
963: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
964: @*/
965: PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
966: {
972: ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);
973: return(0);
974: }
976: /*@
977: DMLabelHasStratum - Determine whether points exist with the given value
979: Input Parameters:
980: + label - the DMLabel
981: - value - the stratum value
983: Output Paramater:
984: . exists - Flag saying whether points exist
986: Level: intermediate
988: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
989: @*/
990: PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists)
991: {
992: PetscInt v;
998: DMLabelLookupStratum(label, value, &v);
999: *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE;
1000: return(0);
1001: }
1003: /*@
1004: DMLabelGetStratumSize - Get the size of a stratum
1006: Input Parameters:
1007: + label - the DMLabel
1008: - value - the stratum value
1010: Output Paramater:
1011: . size - The number of points in the stratum
1013: Level: intermediate
1015: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1016: @*/
1017: PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
1018: {
1019: PetscInt v;
1025: *size = 0;
1026: DMLabelLookupStratum(label, value, &v);
1027: if (v < 0) return(0);
1028: DMLabelMakeValid_Private(label, v);
1029: *size = label->stratumSizes[v];
1030: return(0);
1031: }
1033: /*@
1034: DMLabelGetStratumBounds - Get the largest and smallest point of a stratum
1036: Input Parameters:
1037: + label - the DMLabel
1038: - value - the stratum value
1040: Output Paramaters:
1041: + start - the smallest point in the stratum
1042: - end - the largest point in the stratum
1044: Level: intermediate
1046: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1047: @*/
1048: PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
1049: {
1050: PetscInt v, min, max;
1057: DMLabelLookupStratum(label, value, &v);
1058: if (v < 0) return(0);
1059: DMLabelMakeValid_Private(label, v);
1060: if (label->stratumSizes[v] <= 0) return(0);
1061: ISGetMinMax(label->points[v], &min, &max);
1062: if (start) *start = min;
1063: if (end) *end = max+1;
1064: return(0);
1065: }
1067: /*@
1068: DMLabelGetStratumIS - Get an IS with the stratum points
1070: Input Parameters:
1071: + label - the DMLabel
1072: - value - the stratum value
1074: Output Paramater:
1075: . points - The stratum points
1077: Level: intermediate
1079: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1080: @*/
1081: PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
1082: {
1083: PetscInt v;
1089: *points = NULL;
1090: DMLabelLookupStratum(label, value, &v);
1091: if (v < 0) return(0);
1092: DMLabelMakeValid_Private(label, v);
1093: PetscObjectReference((PetscObject) label->points[v]);
1094: *points = label->points[v];
1095: return(0);
1096: }
1098: /*@
1099: DMLabelSetStratumIS - Set the stratum points using an IS
1101: Input Parameters:
1102: + label - the DMLabel
1103: . value - the stratum value
1104: - points - The stratum points
1106: Level: intermediate
1108: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1109: @*/
1110: PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is)
1111: {
1112: PetscInt v;
1118: DMLabelLookupAddStratum(label, value, &v);
1119: if (is == label->points[v]) return(0);
1120: DMLabelClearStratum(label, value);
1121: ISGetLocalSize(is, &(label->stratumSizes[v]));
1122: PetscObjectReference((PetscObject)is);
1123: ISDestroy(&(label->points[v]));
1124: label->points[v] = is;
1125: label->validIS[v] = PETSC_TRUE;
1126: if (label->bt) {
1127: const PetscInt *points;
1128: PetscInt p;
1130: ISGetIndices(is,&points);
1131: for (p = 0; p < label->stratumSizes[v]; ++p) {
1132: const PetscInt point = points[p];
1134: 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);
1135: PetscBTSet(label->bt, point - label->pStart);
1136: }
1137: }
1138: return(0);
1139: }
1141: /*@
1142: DMLabelClearStratum - Remove a stratum
1144: Input Parameters:
1145: + label - the DMLabel
1146: - value - the stratum value
1148: Level: intermediate
1150: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1151: @*/
1152: PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
1153: {
1154: PetscInt v;
1159: DMLabelLookupStratum(label, value, &v);
1160: if (v < 0) return(0);
1161: if (label->validIS[v]) {
1162: if (label->bt) {
1163: PetscInt i;
1164: const PetscInt *points;
1166: ISGetIndices(label->points[v], &points);
1167: for (i = 0; i < label->stratumSizes[v]; ++i) {
1168: const PetscInt point = points[i];
1170: 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);
1171: PetscBTClear(label->bt, point - label->pStart);
1172: }
1173: ISRestoreIndices(label->points[v], &points);
1174: }
1175: label->stratumSizes[v] = 0;
1176: ISDestroy(&label->points[v]);
1177: ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_OWN_POINTER, &label->points[v]);
1178: PetscObjectSetName((PetscObject) label->points[v], "indices");
1179: } else {
1180: PetscHSetIClear(label->ht[v]);
1181: }
1182: return(0);
1183: }
1185: /*@
1186: DMLabelFilter - Remove all points outside of [start, end)
1188: Input Parameters:
1189: + label - the DMLabel
1190: . start - the first point
1191: - end - the last point
1193: Level: intermediate
1195: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1196: @*/
1197: PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
1198: {
1199: PetscInt v;
1204: DMLabelDestroyIndex(label);
1205: DMLabelMakeAllValid_Private(label);
1206: for (v = 0; v < label->numStrata; ++v) {
1207: PetscInt off, q;
1208: const PetscInt *points;
1209: PetscInt numPointsNew = 0, *pointsNew = NULL;
1211: ISGetIndices(label->points[v], &points);
1212: for (q = 0; q < label->stratumSizes[v]; ++q)
1213: if (points[q] >= start && points[q] < end)
1214: numPointsNew++;
1215: PetscMalloc1(numPointsNew, &pointsNew);
1216: for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) {
1217: if (points[q] >= start && points[q] < end)
1218: pointsNew[off++] = points[q];
1219: }
1220: ISRestoreIndices(label->points[v],&points);
1222: label->stratumSizes[v] = numPointsNew;
1223: ISDestroy(&label->points[v]);
1224: ISCreateGeneral(PETSC_COMM_SELF,numPointsNew, pointsNew, PETSC_OWN_POINTER, &label->points[v]);
1225: PetscObjectSetName((PetscObject) label->points[v], "indices");
1226: }
1227: DMLabelCreateIndex(label, start, end);
1228: return(0);
1229: }
1231: /*@
1232: DMLabelPermute - Create a new label with permuted points
1234: Input Parameters:
1235: + label - the DMLabel
1236: - permutation - the point permutation
1238: Output Parameter:
1239: . labelnew - the new label containing the permuted points
1241: Level: intermediate
1243: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1244: @*/
1245: PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
1246: {
1247: const PetscInt *perm;
1248: PetscInt numValues, numPoints, v, q;
1249: PetscErrorCode ierr;
1254: DMLabelMakeAllValid_Private(label);
1255: DMLabelDuplicate(label, labelNew);
1256: DMLabelGetNumValues(*labelNew, &numValues);
1257: ISGetLocalSize(permutation, &numPoints);
1258: ISGetIndices(permutation, &perm);
1259: for (v = 0; v < numValues; ++v) {
1260: const PetscInt size = (*labelNew)->stratumSizes[v];
1261: const PetscInt *points;
1262: PetscInt *pointsNew;
1264: ISGetIndices((*labelNew)->points[v],&points);
1265: PetscMalloc1(size,&pointsNew);
1266: for (q = 0; q < size; ++q) {
1267: const PetscInt point = points[q];
1269: 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);
1270: pointsNew[q] = perm[point];
1271: }
1272: ISRestoreIndices((*labelNew)->points[v],&points);
1273: PetscSortInt(size, pointsNew);
1274: ISDestroy(&((*labelNew)->points[v]));
1275: if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) {
1276: ISCreateStride(PETSC_COMM_SELF,size,pointsNew[0],1,&((*labelNew)->points[v]));
1277: PetscFree(pointsNew);
1278: } else {
1279: ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));
1280: }
1281: PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");
1282: }
1283: ISRestoreIndices(permutation, &perm);
1284: if (label->bt) {
1285: PetscBTDestroy(&label->bt);
1286: DMLabelCreateIndex(label, label->pStart, label->pEnd);
1287: }
1288: return(0);
1289: }
1291: PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata)
1292: {
1293: MPI_Comm comm;
1294: PetscInt s, l, nroots, nleaves, dof, offset, size;
1295: PetscInt *remoteOffsets, *rootStrata, *rootIdx;
1296: PetscSection rootSection;
1297: PetscSF labelSF;
1301: if (label) {DMLabelMakeAllValid_Private(label);}
1302: PetscObjectGetComm((PetscObject)sf, &comm);
1303: /* Build a section of stratum values per point, generate the according SF
1304: and distribute point-wise stratum values to leaves. */
1305: PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);
1306: PetscSectionCreate(comm, &rootSection);
1307: PetscSectionSetChart(rootSection, 0, nroots);
1308: if (label) {
1309: for (s = 0; s < label->numStrata; ++s) {
1310: const PetscInt *points;
1312: ISGetIndices(label->points[s], &points);
1313: for (l = 0; l < label->stratumSizes[s]; l++) {
1314: PetscSectionGetDof(rootSection, points[l], &dof);
1315: PetscSectionSetDof(rootSection, points[l], dof+1);
1316: }
1317: ISRestoreIndices(label->points[s], &points);
1318: }
1319: }
1320: PetscSectionSetUp(rootSection);
1321: /* Create a point-wise array of stratum values */
1322: PetscSectionGetStorageSize(rootSection, &size);
1323: PetscMalloc1(size, &rootStrata);
1324: PetscCalloc1(nroots, &rootIdx);
1325: if (label) {
1326: for (s = 0; s < label->numStrata; ++s) {
1327: const PetscInt *points;
1329: ISGetIndices(label->points[s], &points);
1330: for (l = 0; l < label->stratumSizes[s]; l++) {
1331: const PetscInt p = points[l];
1332: PetscSectionGetOffset(rootSection, p, &offset);
1333: rootStrata[offset+rootIdx[p]++] = label->stratumValues[s];
1334: }
1335: ISRestoreIndices(label->points[s], &points);
1336: }
1337: }
1338: /* Build SF that maps label points to remote processes */
1339: PetscSectionCreate(comm, leafSection);
1340: PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);
1341: PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);
1342: PetscFree(remoteOffsets);
1343: /* Send the strata for each point over the derived SF */
1344: PetscSectionGetStorageSize(*leafSection, &size);
1345: PetscMalloc1(size, leafStrata);
1346: PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);
1347: PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);
1348: /* Clean up */
1349: PetscFree(rootStrata);
1350: PetscFree(rootIdx);
1351: PetscSectionDestroy(&rootSection);
1352: PetscSFDestroy(&labelSF);
1353: return(0);
1354: }
1356: /*@
1357: DMLabelDistribute - Create a new label pushed forward over the PetscSF
1359: Input Parameters:
1360: + label - the DMLabel
1361: - sf - the map from old to new distribution
1363: Output Parameter:
1364: . labelnew - the new redistributed label
1366: Level: intermediate
1368: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1369: @*/
1370: PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew)
1371: {
1372: MPI_Comm comm;
1373: PetscSection leafSection;
1374: PetscInt p, pStart, pEnd, s, size, dof, offset, stratum;
1375: PetscInt *leafStrata, *strataIdx;
1376: PetscInt **points;
1377: const char *lname = NULL;
1378: char *name;
1379: PetscInt nameSize;
1380: PetscHSetI stratumHash;
1381: size_t len = 0;
1382: PetscMPIInt rank;
1387: if (label) {
1389: DMLabelMakeAllValid_Private(label);
1390: }
1391: PetscObjectGetComm((PetscObject)sf, &comm);
1392: MPI_Comm_rank(comm, &rank);
1393: /* Bcast name */
1394: if (!rank) {
1395: PetscObjectGetName((PetscObject) label, &lname);
1396: PetscStrlen(lname, &len);
1397: }
1398: nameSize = len;
1399: MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);
1400: PetscMalloc1(nameSize+1, &name);
1401: if (!rank) {PetscMemcpy(name, lname, nameSize+1);}
1402: MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);
1403: DMLabelCreate(PETSC_COMM_SELF, name, labelNew);
1404: PetscFree(name);
1405: /* Bcast defaultValue */
1406: if (!rank) (*labelNew)->defaultValue = label->defaultValue;
1407: MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);
1408: /* Distribute stratum values over the SF and get the point mapping on the receiver */
1409: DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);
1410: /* Determine received stratum values and initialise new label*/
1411: PetscHSetICreate(&stratumHash);
1412: PetscSectionGetStorageSize(leafSection, &size);
1413: for (p = 0; p < size; ++p) {PetscHSetIAdd(stratumHash, leafStrata[p]);}
1414: PetscHSetIGetSize(stratumHash, &(*labelNew)->numStrata);
1415: PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);
1416: for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE;
1417: PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);
1418: /* Turn leafStrata into indices rather than stratum values */
1419: offset = 0;
1420: PetscHSetIGetElems(stratumHash, &offset, (*labelNew)->stratumValues);
1421: PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);
1422: for (p = 0; p < size; ++p) {
1423: for (s = 0; s < (*labelNew)->numStrata; ++s) {
1424: if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;}
1425: }
1426: }
1427: /* Rebuild the point strata on the receiver */
1428: PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);
1429: PetscSectionGetChart(leafSection, &pStart, &pEnd);
1430: for (p=pStart; p<pEnd; p++) {
1431: PetscSectionGetDof(leafSection, p, &dof);
1432: PetscSectionGetOffset(leafSection, p, &offset);
1433: for (s=0; s<dof; s++) {
1434: (*labelNew)->stratumSizes[leafStrata[offset+s]]++;
1435: }
1436: }
1437: PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);
1438: PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);
1439: PetscMalloc1((*labelNew)->numStrata,&points);
1440: for (s = 0; s < (*labelNew)->numStrata; ++s) {
1441: PetscHSetICreate(&(*labelNew)->ht[s]);
1442: PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));
1443: }
1444: /* Insert points into new strata */
1445: PetscCalloc1((*labelNew)->numStrata, &strataIdx);
1446: PetscSectionGetChart(leafSection, &pStart, &pEnd);
1447: for (p=pStart; p<pEnd; p++) {
1448: PetscSectionGetDof(leafSection, p, &dof);
1449: PetscSectionGetOffset(leafSection, p, &offset);
1450: for (s=0; s<dof; s++) {
1451: stratum = leafStrata[offset+s];
1452: points[stratum][strataIdx[stratum]++] = p;
1453: }
1454: }
1455: for (s = 0; s < (*labelNew)->numStrata; s++) {
1456: ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));
1457: PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");
1458: }
1459: PetscFree(points);
1460: PetscHSetIDestroy(&stratumHash);
1461: PetscFree(leafStrata);
1462: PetscFree(strataIdx);
1463: PetscSectionDestroy(&leafSection);
1464: return(0);
1465: }
1467: /*@
1468: DMLabelGather - Gather all label values from leafs into roots
1470: Input Parameters:
1471: + label - the DMLabel
1472: - sf - the Star Forest point communication map
1474: Output Parameters:
1475: . labelNew - the new DMLabel with localised leaf values
1477: Level: developer
1479: Note: This is the inverse operation to DMLabelDistribute.
1481: .seealso: DMLabelDistribute()
1482: @*/
1483: PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew)
1484: {
1485: MPI_Comm comm;
1486: PetscSection rootSection;
1487: PetscSF sfLabel;
1488: PetscSFNode *rootPoints, *leafPoints;
1489: PetscInt p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset;
1490: const PetscInt *rootDegree, *ilocal;
1491: PetscInt *rootStrata;
1492: const char *lname;
1493: char *name;
1494: PetscInt nameSize;
1495: size_t len = 0;
1496: PetscMPIInt rank, size;
1502: PetscObjectGetComm((PetscObject)sf, &comm);
1503: MPI_Comm_rank(comm, &rank);
1504: MPI_Comm_size(comm, &size);
1505: /* Bcast name */
1506: if (!rank) {
1507: PetscObjectGetName((PetscObject) label, &lname);
1508: PetscStrlen(lname, &len);
1509: }
1510: nameSize = len;
1511: MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);
1512: PetscMalloc1(nameSize+1, &name);
1513: if (!rank) {PetscMemcpy(name, lname, nameSize+1);}
1514: MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);
1515: DMLabelCreate(PETSC_COMM_SELF, name, labelNew);
1516: PetscFree(name);
1517: /* Gather rank/index pairs of leaves into local roots to build
1518: an inverse, multi-rooted SF. Note that this ignores local leaf
1519: indexing due to the use of the multiSF in PetscSFGather. */
1520: PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);
1521: PetscMalloc1(nroots, &leafPoints);
1522: for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1;
1523: for (p = 0; p < nleaves; p++) {
1524: PetscInt ilp = ilocal ? ilocal[p] : p;
1526: leafPoints[ilp].index = ilp;
1527: leafPoints[ilp].rank = rank;
1528: }
1529: PetscSFComputeDegreeBegin(sf, &rootDegree);
1530: PetscSFComputeDegreeEnd(sf, &rootDegree);
1531: for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p];
1532: PetscMalloc1(nmultiroots, &rootPoints);
1533: PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);
1534: PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);
1535: PetscSFCreate(comm,& sfLabel);
1536: PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);
1537: /* Migrate label over inverted SF to pull stratum values at leaves into roots. */
1538: DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);
1539: /* Rebuild the point strata on the receiver */
1540: for (p = 0, idx = 0; p < nroots; p++) {
1541: for (d = 0; d < rootDegree[p]; d++) {
1542: PetscSectionGetDof(rootSection, idx+d, &dof);
1543: PetscSectionGetOffset(rootSection, idx+d, &offset);
1544: for (s = 0; s < dof; s++) {DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);}
1545: }
1546: idx += rootDegree[p];
1547: }
1548: PetscFree(leafPoints);
1549: PetscFree(rootStrata);
1550: PetscSectionDestroy(&rootSection);
1551: PetscSFDestroy(&sfLabel);
1552: return(0);
1553: }
1555: /*@
1556: DMLabelConvertToSection - Make a PetscSection/IS pair that encodes the label
1558: Input Parameter:
1559: . label - the DMLabel
1561: Output Parameters:
1562: + section - the section giving offsets for each stratum
1563: - is - An IS containing all the label points
1565: Level: developer
1567: .seealso: DMLabelDistribute()
1568: @*/
1569: PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is)
1570: {
1571: IS vIS;
1572: const PetscInt *values;
1573: PetscInt *points;
1574: PetscInt nV, vS = 0, vE = 0, v, N;
1575: PetscErrorCode ierr;
1579: DMLabelGetNumValues(label, &nV);
1580: DMLabelGetValueIS(label, &vIS);
1581: ISGetIndices(vIS, &values);
1582: if (nV) {vS = values[0]; vE = values[0]+1;}
1583: for (v = 1; v < nV; ++v) {
1584: vS = PetscMin(vS, values[v]);
1585: vE = PetscMax(vE, values[v]+1);
1586: }
1587: PetscSectionCreate(PETSC_COMM_SELF, section);
1588: PetscSectionSetChart(*section, vS, vE);
1589: for (v = 0; v < nV; ++v) {
1590: PetscInt n;
1592: DMLabelGetStratumSize(label, values[v], &n);
1593: PetscSectionSetDof(*section, values[v], n);
1594: }
1595: PetscSectionSetUp(*section);
1596: PetscSectionGetStorageSize(*section, &N);
1597: PetscMalloc1(N, &points);
1598: for (v = 0; v < nV; ++v) {
1599: IS is;
1600: const PetscInt *spoints;
1601: PetscInt dof, off, p;
1603: PetscSectionGetDof(*section, values[v], &dof);
1604: PetscSectionGetOffset(*section, values[v], &off);
1605: DMLabelGetStratumIS(label, values[v], &is);
1606: ISGetIndices(is, &spoints);
1607: for (p = 0; p < dof; ++p) points[off+p] = spoints[p];
1608: ISRestoreIndices(is, &spoints);
1609: ISDestroy(&is);
1610: }
1611: ISRestoreIndices(vIS, &values);
1612: ISDestroy(&vIS);
1613: ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);
1614: return(0);
1615: }
1617: /*@
1618: PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using
1619: the local section and an SF describing the section point overlap.
1621: Input Parameters:
1622: + s - The PetscSection for the local field layout
1623: . sf - The SF describing parallel layout of the section points
1624: . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1625: . label - The label specifying the points
1626: - labelValue - The label stratum specifying the points
1628: Output Parameter:
1629: . gsection - The PetscSection for the global field layout
1631: Note: This gives negative sizes and offsets to points not owned by this process
1633: Level: developer
1635: .seealso: PetscSectionCreate()
1636: @*/
1637: PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection)
1638: {
1639: PetscInt *neg = NULL, *tmpOff = NULL;
1640: PetscInt pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots;
1647: PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);
1648: PetscSectionGetChart(s, &pStart, &pEnd);
1649: PetscSectionSetChart(*gsection, pStart, pEnd);
1650: PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);
1651: if (nroots >= 0) {
1652: if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart);
1653: PetscCalloc1(nroots, &neg);
1654: if (nroots > pEnd-pStart) {
1655: PetscCalloc1(nroots, &tmpOff);
1656: } else {
1657: tmpOff = &(*gsection)->atlasDof[-pStart];
1658: }
1659: }
1660: /* Mark ghost points with negative dof */
1661: for (p = pStart; p < pEnd; ++p) {
1662: PetscInt value;
1664: DMLabelGetValue(label, p, &value);
1665: if (value != labelValue) continue;
1666: PetscSectionGetDof(s, p, &dof);
1667: PetscSectionSetDof(*gsection, p, dof);
1668: PetscSectionGetConstraintDof(s, p, &cdof);
1669: if (!includeConstraints && cdof > 0) {PetscSectionSetConstraintDof(*gsection, p, cdof);}
1670: if (neg) neg[p] = -(dof+1);
1671: }
1672: PetscSectionSetUpBC(*gsection);
1673: if (nroots >= 0) {
1674: PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);
1675: PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);
1676: if (nroots > pEnd-pStart) {
1677: for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
1678: }
1679: }
1680: /* Calculate new sizes, get proccess offset, and calculate point offsets */
1681: for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1682: cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
1683: (*gsection)->atlasOff[p] = off;
1684: off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
1685: }
1686: MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));
1687: globalOff -= off;
1688: for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1689: (*gsection)->atlasOff[p] += globalOff;
1690: if (neg) neg[p] = -((*gsection)->atlasOff[p]+1);
1691: }
1692: /* Put in negative offsets for ghost points */
1693: if (nroots >= 0) {
1694: PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);
1695: PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);
1696: if (nroots > pEnd-pStart) {
1697: for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1698: }
1699: }
1700: if (nroots >= 0 && nroots > pEnd-pStart) {PetscFree(tmpOff);}
1701: PetscFree(neg);
1702: return(0);
1703: }
1705: typedef struct _n_PetscSectionSym_Label
1706: {
1707: DMLabel label;
1708: PetscCopyMode *modes;
1709: PetscInt *sizes;
1710: const PetscInt ***perms;
1711: const PetscScalar ***rots;
1712: PetscInt (*minMaxOrients)[2];
1713: PetscInt numStrata; /* numStrata is only increasing, functions as a state */
1714: } PetscSectionSym_Label;
1716: static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym)
1717: {
1718: PetscInt i, j;
1719: PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
1720: PetscErrorCode ierr;
1723: for (i = 0; i <= sl->numStrata; i++) {
1724: if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) {
1725: for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
1726: if (sl->perms[i]) {PetscFree(sl->perms[i][j]);}
1727: if (sl->rots[i]) {PetscFree(sl->rots[i][j]);}
1728: }
1729: if (sl->perms[i]) {
1730: const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]];
1732: PetscFree(perms);
1733: }
1734: if (sl->rots[i]) {
1735: const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]];
1737: PetscFree(rots);
1738: }
1739: }
1740: }
1741: PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);
1742: DMLabelDestroy(&sl->label);
1743: sl->numStrata = 0;
1744: return(0);
1745: }
1747: static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym)
1748: {
1752: PetscSectionSymLabelReset(sym);
1753: PetscFree(sym->data);
1754: return(0);
1755: }
1757: static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer)
1758: {
1759: PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
1760: PetscBool isAscii;
1761: DMLabel label = sl->label;
1762: const char *name;
1763: PetscErrorCode ierr;
1766: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);
1767: if (isAscii) {
1768: PetscInt i, j, k;
1769: PetscViewerFormat format;
1771: PetscViewerGetFormat(viewer,&format);
1772: if (label) {
1773: PetscViewerGetFormat(viewer,&format);
1774: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1775: PetscViewerASCIIPushTab(viewer);
1776: DMLabelView(label, viewer);
1777: PetscViewerASCIIPopTab(viewer);
1778: } else {
1779: PetscObjectGetName((PetscObject) sl->label, &name);
1780: PetscViewerASCIIPrintf(viewer," Label '%s'\n",name);
1781: }
1782: } else {
1783: PetscViewerASCIIPrintf(viewer, "No label given\n");
1784: }
1785: PetscViewerASCIIPushTab(viewer);
1786: for (i = 0; i <= sl->numStrata; i++) {
1787: PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue;
1789: if (!(sl->perms[i] || sl->rots[i])) {
1790: PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);
1791: } else {
1792: PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);
1793: PetscViewerASCIIPushTab(viewer);
1794: PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);
1795: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1796: PetscViewerASCIIPushTab(viewer);
1797: for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
1798: if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) {
1799: PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);
1800: } else {
1801: PetscInt tab;
1803: PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);
1804: PetscViewerASCIIPushTab(viewer);
1805: PetscViewerASCIIGetTab(viewer,&tab);
1806: if (sl->perms[i] && sl->perms[i][j]) {
1807: PetscViewerASCIIPrintf(viewer,"Permutation:");
1808: PetscViewerASCIISetTab(viewer,0);
1809: for (k = 0; k < sl->sizes[i]; k++) {PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);}
1810: PetscViewerASCIIPrintf(viewer,"\n");
1811: PetscViewerASCIISetTab(viewer,tab);
1812: }
1813: if (sl->rots[i] && sl->rots[i][j]) {
1814: PetscViewerASCIIPrintf(viewer,"Rotations: ");
1815: PetscViewerASCIISetTab(viewer,0);
1816: #if defined(PETSC_USE_COMPLEX)
1817: for (k = 0; k < sl->sizes[i]; k++) {PetscViewerASCIIPrintf(viewer," %+f+i*%+f",PetscRealPart(sl->rots[i][j][k]),PetscImaginaryPart(sl->rots[i][j][k]));}
1818: #else
1819: for (k = 0; k < sl->sizes[i]; k++) {PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);}
1820: #endif
1821: PetscViewerASCIIPrintf(viewer,"\n");
1822: PetscViewerASCIISetTab(viewer,tab);
1823: }
1824: PetscViewerASCIIPopTab(viewer);
1825: }
1826: }
1827: PetscViewerASCIIPopTab(viewer);
1828: }
1829: PetscViewerASCIIPopTab(viewer);
1830: }
1831: }
1832: PetscViewerASCIIPopTab(viewer);
1833: }
1834: return(0);
1835: }
1837: /*@
1838: PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries
1840: Logically collective on sym
1842: Input parameters:
1843: + sym - the section symmetries
1844: - label - the DMLabel describing the types of points
1846: Level: developer:
1848: .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms()
1849: @*/
1850: PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label)
1851: {
1852: PetscSectionSym_Label *sl;
1853: PetscErrorCode ierr;
1857: sl = (PetscSectionSym_Label *) sym->data;
1858: if (sl->label && sl->label != label) {PetscSectionSymLabelReset(sym);}
1859: if (label) {
1860: sl->label = label;
1861: PetscObjectReference((PetscObject) label);
1862: DMLabelGetNumValues(label,&sl->numStrata);
1863: PetscMalloc5(sl->numStrata+1,&sl->modes,sl->numStrata+1,&sl->sizes,sl->numStrata+1,&sl->perms,sl->numStrata+1,&sl->rots,sl->numStrata+1,&sl->minMaxOrients);
1864: PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));
1865: PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));
1866: PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));
1867: PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));
1868: PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));
1869: }
1870: return(0);
1871: }
1873: /*@C
1874: PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum
1876: Logically collective on PetscSectionSym
1878: InputParameters:
1879: + sys - the section symmetries
1880: . stratum - the stratum value in the label that we are assigning symmetries for
1881: . size - the number of dofs for points in the stratum of the label
1882: . minOrient - the smallest orientation for a point in this stratum
1883: . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient))
1884: . mode - how sym should copy the perms and rots arrays
1885: . perms - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation. A NULL permutation is the identity
1886: + rots - NULL if there are no rotations, or (maxOrient - minOrient) sets of rotations, one for each orientation. A NULL set of orientations is the identity
1888: Level: developer
1890: .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel()
1891: @*/
1892: PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots)
1893: {
1894: PetscSectionSym_Label *sl;
1895: const char *name;
1896: PetscInt i, j, k;
1897: PetscErrorCode ierr;
1901: sl = (PetscSectionSym_Label *) sym->data;
1902: if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet");
1903: for (i = 0; i <= sl->numStrata; i++) {
1904: PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue;
1906: if (stratum == value) break;
1907: }
1908: PetscObjectGetName((PetscObject) sl->label, &name);
1909: if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,name);
1910: sl->sizes[i] = size;
1911: sl->modes[i] = mode;
1912: sl->minMaxOrients[i][0] = minOrient;
1913: sl->minMaxOrients[i][1] = maxOrient;
1914: if (mode == PETSC_COPY_VALUES) {
1915: if (perms) {
1916: PetscInt **ownPerms;
1918: PetscCalloc1(maxOrient - minOrient,&ownPerms);
1919: for (j = 0; j < maxOrient-minOrient; j++) {
1920: if (perms[j]) {
1921: PetscMalloc1(size,&ownPerms[j]);
1922: for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];}
1923: }
1924: }
1925: sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient];
1926: }
1927: if (rots) {
1928: PetscScalar **ownRots;
1930: PetscCalloc1(maxOrient - minOrient,&ownRots);
1931: for (j = 0; j < maxOrient-minOrient; j++) {
1932: if (rots[j]) {
1933: PetscMalloc1(size,&ownRots[j]);
1934: for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];}
1935: }
1936: }
1937: sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient];
1938: }
1939: } else {
1940: sl->perms[i] = perms ? &perms[-minOrient] : NULL;
1941: sl->rots[i] = rots ? &rots[-minOrient] : NULL;
1942: }
1943: return(0);
1944: }
1946: static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots)
1947: {
1948: PetscInt i, j, numStrata;
1949: PetscSectionSym_Label *sl;
1950: DMLabel label;
1951: PetscErrorCode ierr;
1954: sl = (PetscSectionSym_Label *) sym->data;
1955: numStrata = sl->numStrata;
1956: label = sl->label;
1957: for (i = 0; i < numPoints; i++) {
1958: PetscInt point = points[2*i];
1959: PetscInt ornt = points[2*i+1];
1961: for (j = 0; j < numStrata; j++) {
1962: if (label->validIS[j]) {
1963: PetscInt k;
1965: ISLocate(label->points[j],point,&k);
1966: if (k >= 0) break;
1967: } else {
1968: PetscBool has;
1970: PetscHSetIHas(label->ht[j], point, &has);
1971: if (has) break;
1972: }
1973: }
1974: if ((sl->minMaxOrients[j][1] > sl->minMaxOrients[j][0]) && (ornt < sl->minMaxOrients[j][0] || ornt >= sl->minMaxOrients[j][1])) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"point %D orientation %D not in range [%D, %D) for stratum %D",point,ornt,sl->minMaxOrients[j][0],sl->minMaxOrients[j][1],j < numStrata ? label->stratumValues[j] : label->defaultValue);
1975: if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;}
1976: if (rots) {rots[i] = sl->rots[j] ? sl->rots[j][ornt] : NULL;}
1977: }
1978: return(0);
1979: }
1981: PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym)
1982: {
1983: PetscSectionSym_Label *sl;
1984: PetscErrorCode ierr;
1987: PetscNewLog(sym,&sl);
1988: sym->ops->getpoints = PetscSectionSymGetPoints_Label;
1989: sym->ops->view = PetscSectionSymView_Label;
1990: sym->ops->destroy = PetscSectionSymDestroy_Label;
1991: sym->data = (void *) sl;
1992: return(0);
1993: }
1995: /*@
1996: PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label
1998: Collective
2000: Input Parameters:
2001: + comm - the MPI communicator for the new symmetry
2002: - label - the label defining the strata
2004: Output Parameters:
2005: . sym - the section symmetries
2007: Level: developer
2009: .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()
2010: @*/
2011: PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym)
2012: {
2013: PetscErrorCode ierr;
2016: DMInitializePackage();
2017: PetscSectionSymCreate(comm,sym);
2018: PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);
2019: PetscSectionSymLabelSetLabel(*sym,label);
2020: return(0);
2021: }