Actual source code: dmlabel.c
1: #include <petscdm.h>
2: #include <petsc/private/dmlabelimpl.h>
3: #include <petsc/private/sectionimpl.h>
4: #include <petscsf.h>
5: #include <petscsection.h>
7: /*@C
8: DMLabelCreate - Create a DMLabel object, which is a multimap
10: Collective
12: Input parameters:
13: + comm - The communicator, usually PETSC_COMM_SELF
14: - name - The label name
16: Output parameter:
17: . label - The DMLabel
19: Level: beginner
21: Notes:
22: The label name is actually usual PetscObject name.
23: One can get/set it with PetscObjectGetName()/PetscObjectSetName().
25: .seealso: DMLabelDestroy()
26: @*/
27: PetscErrorCode DMLabelCreate(MPI_Comm comm, const char name[], DMLabel *label)
28: {
33: DMInitializePackage();
35: PetscHeaderCreate(*label,DMLABEL_CLASSID,"DMLabel","DMLabel","DM",comm,DMLabelDestroy,DMLabelView);
37: (*label)->numStrata = 0;
38: (*label)->defaultValue = -1;
39: (*label)->stratumValues = NULL;
40: (*label)->validIS = NULL;
41: (*label)->stratumSizes = NULL;
42: (*label)->points = NULL;
43: (*label)->ht = NULL;
44: (*label)->pStart = -1;
45: (*label)->pEnd = -1;
46: (*label)->bt = NULL;
47: PetscHMapICreate(&(*label)->hmap);
48: PetscObjectSetName((PetscObject) *label, name);
49: return(0);
50: }
52: /*
53: DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format
55: Not collective
57: Input parameter:
58: + label - The DMLabel
59: - v - The stratum value
61: Output parameter:
62: . label - The DMLabel with stratum in sorted list format
64: Level: developer
66: .seealso: DMLabelCreate()
67: */
68: static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v)
69: {
70: IS is;
71: PetscInt off = 0, *pointArray, p;
74: if (PetscLikely(v >= 0 && v < label->numStrata) && label->validIS[v]) return 0;
76: 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);
77: PetscHSetIGetSize(label->ht[v], &label->stratumSizes[v]);
78: PetscMalloc1(label->stratumSizes[v], &pointArray);
79: PetscHSetIGetElems(label->ht[v], &off, pointArray);
80: PetscHSetIClear(label->ht[v]);
81: PetscSortInt(label->stratumSizes[v], pointArray);
82: if (label->bt) {
83: for (p = 0; p < label->stratumSizes[v]; ++p) {
84: const PetscInt point = pointArray[p];
85: 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);
86: PetscBTSet(label->bt, point - label->pStart);
87: }
88: }
89: if (label->stratumSizes[v] > 0 && pointArray[label->stratumSizes[v]-1] == pointArray[0] + label->stratumSizes[v]-1) {
90: ISCreateStride(PETSC_COMM_SELF, label->stratumSizes[v], pointArray[0], 1, &is);
91: PetscFree(pointArray);
92: } else {
93: ISCreateGeneral(PETSC_COMM_SELF, label->stratumSizes[v], pointArray, PETSC_OWN_POINTER, &is);
94: }
95: PetscObjectSetName((PetscObject) is, "indices");
96: label->points[v] = is;
97: label->validIS[v] = PETSC_TRUE;
98: PetscObjectStateIncrease((PetscObject) label);
99: return(0);
100: }
102: /*
103: DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format
105: Not collective
107: Input parameter:
108: . label - The DMLabel
110: Output parameter:
111: . label - The DMLabel with all strata in sorted list format
113: Level: developer
115: .seealso: DMLabelCreate()
116: */
117: static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label)
118: {
119: PetscInt v;
123: for (v = 0; v < label->numStrata; v++) {
124: DMLabelMakeValid_Private(label, v);
125: }
126: return(0);
127: }
129: /*
130: DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format
132: Not collective
134: Input parameter:
135: + label - The DMLabel
136: - v - The stratum value
138: Output parameter:
139: . label - The DMLabel with stratum in hash format
141: Level: developer
143: .seealso: DMLabelCreate()
144: */
145: static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v)
146: {
147: PetscInt p;
148: const PetscInt *points;
151: if (PetscLikely(v >= 0 && v < label->numStrata) && !label->validIS[v]) return 0;
153: 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);
154: if (label->points[v]) {
155: ISGetIndices(label->points[v], &points);
156: for (p = 0; p < label->stratumSizes[v]; ++p) {
157: PetscHSetIAdd(label->ht[v], points[p]);
158: }
159: ISRestoreIndices(label->points[v],&points);
160: ISDestroy(&(label->points[v]));
161: }
162: label->validIS[v] = PETSC_FALSE;
163: return(0);
164: }
166: #if !defined(DMLABEL_LOOKUP_THRESHOLD)
167: #define DMLABEL_LOOKUP_THRESHOLD 16
168: #endif
170: PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupStratum(DMLabel label, PetscInt value, PetscInt *index)
171: {
172: PetscInt v;
176: *index = -1;
177: if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD) {
178: for (v = 0; v < label->numStrata; ++v)
179: if (label->stratumValues[v] == value) {*index = v; break;}
180: } else {
181: PetscHMapIGet(label->hmap, value, index);
182: }
183: if (PetscDefined(USE_DEBUG)) { /* Check strata hash map consistency */
184: PetscInt len, loc = -1;
185: PetscHMapIGetSize(label->hmap, &len);
186: if (len != label->numStrata) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent strata hash map size");
187: if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD) {
188: PetscHMapIGet(label->hmap, value, &loc);
189: } else {
190: for (v = 0; v < label->numStrata; ++v)
191: if (label->stratumValues[v] == value) {loc = v; break;}
192: }
193: if (loc != *index) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent strata hash map lookup");
194: }
195: return(0);
196: }
198: PETSC_STATIC_INLINE PetscErrorCode DMLabelNewStratum(DMLabel label, PetscInt value, PetscInt *index)
199: {
200: PetscInt v;
201: PetscInt *tmpV;
202: PetscInt *tmpS;
203: PetscHSetI *tmpH, ht;
204: IS *tmpP, is;
205: PetscBool *tmpB;
206: PetscHMapI hmap = label->hmap;
210: v = label->numStrata;
211: tmpV = label->stratumValues;
212: tmpS = label->stratumSizes;
213: tmpH = label->ht;
214: tmpP = label->points;
215: tmpB = label->validIS;
216: { /* TODO: PetscRealloc() is broken, use malloc+memcpy+free */
217: PetscInt *oldV = tmpV;
218: PetscInt *oldS = tmpS;
219: PetscHSetI *oldH = tmpH;
220: IS *oldP = tmpP;
221: PetscBool *oldB = tmpB;
222: PetscMalloc((v+1)*sizeof(*tmpV), &tmpV);
223: PetscMalloc((v+1)*sizeof(*tmpS), &tmpS);
224: PetscMalloc((v+1)*sizeof(*tmpH), &tmpH);
225: PetscMalloc((v+1)*sizeof(*tmpP), &tmpP);
226: PetscMalloc((v+1)*sizeof(*tmpB), &tmpB);
227: PetscArraycpy(tmpV, oldV, v);
228: PetscArraycpy(tmpS, oldS, v);
229: PetscArraycpy(tmpH, oldH, v);
230: PetscArraycpy(tmpP, oldP, v);
231: PetscArraycpy(tmpB, oldB, v);
232: PetscFree(oldV);
233: PetscFree(oldS);
234: PetscFree(oldH);
235: PetscFree(oldP);
236: PetscFree(oldB);
237: }
238: label->numStrata = v+1;
239: label->stratumValues = tmpV;
240: label->stratumSizes = tmpS;
241: label->ht = tmpH;
242: label->points = tmpP;
243: label->validIS = tmpB;
244: PetscHSetICreate(&ht);
245: ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);
246: PetscHMapISet(hmap, value, v);
247: tmpV[v] = value;
248: tmpS[v] = 0;
249: tmpH[v] = ht;
250: tmpP[v] = is;
251: tmpB[v] = PETSC_TRUE;
252: PetscObjectStateIncrease((PetscObject) label);
253: *index = v;
254: return(0);
255: }
257: PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupAddStratum(DMLabel label, PetscInt value, PetscInt *index)
258: {
261: DMLabelLookupStratum(label, value, index);
262: if (*index < 0) {DMLabelNewStratum(label, value, index);}
263: return(0);
264: }
266: /*@
267: DMLabelAddStratum - Adds a new stratum value in a DMLabel
269: Input Parameters:
270: + label - The DMLabel
271: - value - The stratum value
273: Level: beginner
275: .seealso: DMLabelCreate(), DMLabelDestroy()
276: @*/
277: PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value)
278: {
279: PetscInt v;
284: DMLabelLookupAddStratum(label, value, &v);
285: return(0);
286: }
288: /*@
289: DMLabelAddStrata - Adds new stratum values in a DMLabel
291: Not collective
293: Input Parameters:
294: + label - The DMLabel
295: . numStrata - The number of stratum values
296: - stratumValues - The stratum values
298: Level: beginner
300: .seealso: DMLabelCreate(), DMLabelDestroy()
301: @*/
302: PetscErrorCode DMLabelAddStrata(DMLabel label, PetscInt numStrata, const PetscInt stratumValues[])
303: {
304: PetscInt *values, v;
310: PetscMalloc1(numStrata, &values);
311: PetscArraycpy(values, stratumValues, numStrata);
312: PetscSortRemoveDupsInt(&numStrata, values);
313: if (!label->numStrata) { /* Fast preallocation */
314: PetscInt *tmpV;
315: PetscInt *tmpS;
316: PetscHSetI *tmpH, ht;
317: IS *tmpP, is;
318: PetscBool *tmpB;
319: PetscHMapI hmap = label->hmap;
321: PetscMalloc1(numStrata, &tmpV);
322: PetscMalloc1(numStrata, &tmpS);
323: PetscMalloc1(numStrata, &tmpH);
324: PetscMalloc1(numStrata, &tmpP);
325: PetscMalloc1(numStrata, &tmpB);
326: label->numStrata = numStrata;
327: label->stratumValues = tmpV;
328: label->stratumSizes = tmpS;
329: label->ht = tmpH;
330: label->points = tmpP;
331: label->validIS = tmpB;
332: for (v = 0; v < numStrata; ++v) {
333: PetscHSetICreate(&ht);
334: ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);
335: PetscHMapISet(hmap, values[v], v);
336: tmpV[v] = values[v];
337: tmpS[v] = 0;
338: tmpH[v] = ht;
339: tmpP[v] = is;
340: tmpB[v] = PETSC_TRUE;
341: }
342: PetscObjectStateIncrease((PetscObject) label);
343: } else {
344: for (v = 0; v < numStrata; ++v) {
345: DMLabelAddStratum(label, values[v]);
346: }
347: }
348: PetscFree(values);
349: return(0);
350: }
352: /*@
353: DMLabelAddStrataIS - Adds new stratum values in a DMLabel
355: Not collective
357: Input Parameters:
358: + label - The DMLabel
359: - valueIS - Index set with stratum values
361: Level: beginner
363: .seealso: DMLabelCreate(), DMLabelDestroy()
364: @*/
365: PetscErrorCode DMLabelAddStrataIS(DMLabel label, IS valueIS)
366: {
367: PetscInt numStrata;
368: const PetscInt *stratumValues;
374: ISGetLocalSize(valueIS, &numStrata);
375: ISGetIndices(valueIS, &stratumValues);
376: DMLabelAddStrata(label, numStrata, stratumValues);
377: return(0);
378: }
380: static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer)
381: {
382: PetscInt v;
383: PetscMPIInt rank;
387: MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);
388: PetscViewerASCIIPushSynchronized(viewer);
389: if (label) {
390: const char *name;
392: PetscObjectGetName((PetscObject) label, &name);
393: PetscViewerASCIIPrintf(viewer, "Label '%s':\n", name);
394: if (label->bt) {PetscViewerASCIIPrintf(viewer, " Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);}
395: for (v = 0; v < label->numStrata; ++v) {
396: const PetscInt value = label->stratumValues[v];
397: const PetscInt *points;
398: PetscInt p;
400: ISGetIndices(label->points[v], &points);
401: for (p = 0; p < label->stratumSizes[v]; ++p) {
402: PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);
403: }
404: ISRestoreIndices(label->points[v],&points);
405: }
406: }
407: PetscViewerFlush(viewer);
408: PetscViewerASCIIPopSynchronized(viewer);
409: return(0);
410: }
412: /*@C
413: DMLabelView - View the label
415: Collective on viewer
417: Input Parameters:
418: + label - The DMLabel
419: - viewer - The PetscViewer
421: Level: intermediate
423: .seealso: DMLabelCreate(), DMLabelDestroy()
424: @*/
425: PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
426: {
427: PetscBool iascii;
432: if (!viewer) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)label), &viewer);}
434: if (label) {DMLabelMakeAllValid_Private(label);}
435: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
436: if (iascii) {
437: DMLabelView_Ascii(label, viewer);
438: }
439: return(0);
440: }
442: /*@
443: DMLabelReset - Destroys internal data structures in a DMLabel
445: Not collective
447: Input Parameter:
448: . label - The DMLabel
450: Level: beginner
452: .seealso: DMLabelDestroy(), DMLabelCreate()
453: @*/
454: PetscErrorCode DMLabelReset(DMLabel label)
455: {
456: PetscInt v;
461: for (v = 0; v < label->numStrata; ++v) {
462: PetscHSetIDestroy(&label->ht[v]);
463: ISDestroy(&label->points[v]);
464: }
465: label->numStrata = 0;
466: PetscFree(label->stratumValues);
467: PetscFree(label->stratumSizes);
468: PetscFree(label->ht);
469: PetscFree(label->points);
470: PetscFree(label->validIS);
471: label->stratumValues = NULL;
472: label->stratumSizes = NULL;
473: label->ht = NULL;
474: label->points = NULL;
475: label->validIS = NULL;
476: PetscHMapIReset(label->hmap);
477: label->pStart = -1;
478: label->pEnd = -1;
479: PetscBTDestroy(&label->bt);
480: return(0);
481: }
483: /*@
484: DMLabelDestroy - Destroys a DMLabel
486: Collective on label
488: Input Parameter:
489: . label - The DMLabel
491: Level: beginner
493: .seealso: DMLabelReset(), DMLabelCreate()
494: @*/
495: PetscErrorCode DMLabelDestroy(DMLabel *label)
496: {
500: if (!*label) return(0);
502: if (--((PetscObject)(*label))->refct > 0) {*label = NULL; return(0);}
503: DMLabelReset(*label);
504: PetscHMapIDestroy(&(*label)->hmap);
505: PetscHeaderDestroy(label);
506: return(0);
507: }
509: /*@
510: DMLabelDuplicate - Duplicates a DMLabel
512: Collective on label
514: Input Parameter:
515: . label - The DMLabel
517: Output Parameter:
518: . labelnew - location to put new vector
520: Level: intermediate
522: .seealso: DMLabelCreate(), DMLabelDestroy()
523: @*/
524: PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
525: {
526: const char *name;
527: PetscInt v;
532: DMLabelMakeAllValid_Private(label);
533: PetscObjectGetName((PetscObject) label, &name);
534: DMLabelCreate(PetscObjectComm((PetscObject) label), name, labelnew);
536: (*labelnew)->numStrata = label->numStrata;
537: (*labelnew)->defaultValue = label->defaultValue;
538: PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);
539: PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);
540: PetscMalloc1(label->numStrata, &(*labelnew)->ht);
541: PetscMalloc1(label->numStrata, &(*labelnew)->points);
542: PetscMalloc1(label->numStrata, &(*labelnew)->validIS);
543: for (v = 0; v < label->numStrata; ++v) {
544: PetscHSetICreate(&(*labelnew)->ht[v]);
545: (*labelnew)->stratumValues[v] = label->stratumValues[v];
546: (*labelnew)->stratumSizes[v] = label->stratumSizes[v];
547: PetscObjectReference((PetscObject) (label->points[v]));
548: (*labelnew)->points[v] = label->points[v];
549: (*labelnew)->validIS[v] = PETSC_TRUE;
550: }
551: PetscHMapIDestroy(&(*labelnew)->hmap);
552: PetscHMapIDuplicate(label->hmap,&(*labelnew)->hmap);
553: (*labelnew)->pStart = -1;
554: (*labelnew)->pEnd = -1;
555: (*labelnew)->bt = NULL;
556: return(0);
557: }
559: /*@
560: DMLabelComputeIndex - Create an index structure for membership determination, automatically determining the bounds
562: Not collective
564: Input Parameter:
565: . label - The DMLabel
567: Level: intermediate
569: .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
570: @*/
571: PetscErrorCode DMLabelComputeIndex(DMLabel label)
572: {
573: PetscInt pStart = PETSC_MAX_INT, pEnd = -1, v;
578: DMLabelMakeAllValid_Private(label);
579: for (v = 0; v < label->numStrata; ++v) {
580: const PetscInt *points;
581: PetscInt i;
583: ISGetIndices(label->points[v], &points);
584: for (i = 0; i < label->stratumSizes[v]; ++i) {
585: const PetscInt point = points[i];
587: pStart = PetscMin(point, pStart);
588: pEnd = PetscMax(point+1, pEnd);
589: }
590: ISRestoreIndices(label->points[v], &points);
591: }
592: label->pStart = pStart == PETSC_MAX_INT ? -1 : pStart;
593: label->pEnd = pEnd;
594: DMLabelCreateIndex(label, label->pStart, label->pEnd);
595: return(0);
596: }
598: /*@
599: DMLabelCreateIndex - Create an index structure for membership determination
601: Not collective
603: Input Parameters:
604: + label - The DMLabel
605: . pStart - The smallest point
606: - pEnd - The largest point + 1
608: Level: intermediate
610: .seealso: DMLabelHasPoint(), DMLabelComputeIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
611: @*/
612: PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
613: {
614: PetscInt v;
619: DMLabelDestroyIndex(label);
620: DMLabelMakeAllValid_Private(label);
621: label->pStart = pStart;
622: label->pEnd = pEnd;
623: /* This can be hooked into SetValue(), ClearValue(), etc. for updating */
624: PetscBTCreate(pEnd - pStart, &label->bt);
625: for (v = 0; v < label->numStrata; ++v) {
626: const PetscInt *points;
627: PetscInt i;
629: ISGetIndices(label->points[v], &points);
630: for (i = 0; i < label->stratumSizes[v]; ++i) {
631: const PetscInt point = points[i];
633: if ((point < pStart) || (point >= pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, pStart, pEnd);
634: PetscBTSet(label->bt, point - pStart);
635: }
636: ISRestoreIndices(label->points[v], &points);
637: }
638: return(0);
639: }
641: /*@
642: DMLabelDestroyIndex - Destroy the index structure
644: Not collective
646: Input Parameter:
647: . label - the DMLabel
649: Level: intermediate
651: .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
652: @*/
653: PetscErrorCode DMLabelDestroyIndex(DMLabel label)
654: {
659: label->pStart = -1;
660: label->pEnd = -1;
661: PetscBTDestroy(&label->bt);
662: return(0);
663: }
665: /*@
666: DMLabelGetBounds - Return the smallest and largest point in the label
668: Not collective
670: Input Parameter:
671: . label - the DMLabel
673: Output Parameters:
674: + pStart - The smallest point
675: - pEnd - The largest point + 1
677: Note: This will compute an index for the label if one does not exist.
679: Level: intermediate
681: .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
682: @*/
683: PetscErrorCode DMLabelGetBounds(DMLabel label, PetscInt *pStart, PetscInt *pEnd)
684: {
689: if ((label->pStart == -1) && (label->pEnd == -1)) {DMLabelComputeIndex(label);}
690: if (pStart) {
692: *pStart = label->pStart;
693: }
694: if (pEnd) {
696: *pEnd = label->pEnd;
697: }
698: return(0);
699: }
701: /*@
702: DMLabelHasValue - Determine whether a label assigns the value to any point
704: Not collective
706: Input Parameters:
707: + label - the DMLabel
708: - value - the value
710: Output Parameter:
711: . contains - Flag indicating whether the label maps this value to any point
713: Level: developer
715: .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue()
716: @*/
717: PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains)
718: {
719: PetscInt v;
725: DMLabelLookupStratum(label, value, &v);
726: *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE;
727: return(0);
728: }
730: /*@
731: DMLabelHasPoint - Determine whether a label assigns a value to a point
733: Not collective
735: Input Parameters:
736: + label - the DMLabel
737: - point - the point
739: Output Parameter:
740: . contains - Flag indicating whether the label maps this point to a value
742: Note: The user must call DMLabelCreateIndex() before this function.
744: Level: developer
746: .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
747: @*/
748: PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
749: {
755: DMLabelMakeAllValid_Private(label);
756: if (PetscDefined(USE_DEBUG)) {
757: if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
758: 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);
759: }
760: *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
761: return(0);
762: }
764: /*@
765: DMLabelStratumHasPoint - Return true if the stratum contains a point
767: Not collective
769: Input Parameters:
770: + label - the DMLabel
771: . value - the stratum value
772: - point - the point
774: Output Parameter:
775: . contains - true if the stratum contains the point
777: Level: intermediate
779: .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue()
780: @*/
781: PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains)
782: {
783: PetscInt v;
789: *contains = PETSC_FALSE;
790: DMLabelLookupStratum(label, value, &v);
791: if (v < 0) return(0);
793: if (label->validIS[v]) {
794: PetscInt i;
796: ISLocate(label->points[v], point, &i);
797: if (i >= 0) *contains = PETSC_TRUE;
798: } else {
799: PetscBool has;
801: PetscHSetIHas(label->ht[v], point, &has);
802: if (has) *contains = PETSC_TRUE;
803: }
804: return(0);
805: }
807: /*@
808: DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
809: When a label is created, it is initialized to -1.
811: Not collective
813: Input parameter:
814: . label - a DMLabel object
816: Output parameter:
817: . defaultValue - the default value
819: Level: beginner
821: .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
822: @*/
823: PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue)
824: {
827: *defaultValue = label->defaultValue;
828: return(0);
829: }
831: /*@
832: DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
833: When a label is created, it is initialized to -1.
835: Not collective
837: Input parameter:
838: . label - a DMLabel object
840: Output parameter:
841: . defaultValue - the default value
843: Level: beginner
845: .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
846: @*/
847: PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue)
848: {
851: label->defaultValue = defaultValue;
852: return(0);
853: }
855: /*@
856: 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())
858: Not collective
860: Input Parameters:
861: + label - the DMLabel
862: - point - the point
864: Output Parameter:
865: . value - The point value, or the default value (-1 by default)
867: Level: intermediate
869: .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
870: @*/
871: PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
872: {
873: PetscInt v;
879: *value = label->defaultValue;
880: for (v = 0; v < label->numStrata; ++v) {
881: if (label->validIS[v]) {
882: PetscInt i;
884: ISLocate(label->points[v], point, &i);
885: if (i >= 0) {
886: *value = label->stratumValues[v];
887: break;
888: }
889: } else {
890: PetscBool has;
892: PetscHSetIHas(label->ht[v], point, &has);
893: if (has) {
894: *value = label->stratumValues[v];
895: break;
896: }
897: }
898: }
899: return(0);
900: }
902: /*@
903: 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.
905: Not collective
907: Input Parameters:
908: + label - the DMLabel
909: . point - the point
910: - value - The point value
912: Level: intermediate
914: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
915: @*/
916: PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
917: {
918: PetscInt v;
923: /* Find label value, add new entry if needed */
924: if (value == label->defaultValue) return(0);
925: DMLabelLookupAddStratum(label, value, &v);
926: /* Set key */
927: DMLabelMakeInvalid_Private(label, v);
928: PetscHSetIAdd(label->ht[v], point);
929: return(0);
930: }
932: /*@
933: DMLabelClearValue - Clear the value a label assigns to a point
935: Not collective
937: Input Parameters:
938: + label - the DMLabel
939: . point - the point
940: - value - The point value
942: Level: intermediate
944: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue()
945: @*/
946: PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
947: {
948: PetscInt v;
953: /* Find label value */
954: DMLabelLookupStratum(label, value, &v);
955: if (v < 0) return(0);
957: if (label->bt) {
958: 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);
959: PetscBTClear(label->bt, point - label->pStart);
960: }
962: /* Delete key */
963: DMLabelMakeInvalid_Private(label, v);
964: PetscHSetIDel(label->ht[v], point);
965: return(0);
966: }
968: /*@
969: DMLabelInsertIS - Set all points in the IS to a value
971: Not collective
973: Input Parameters:
974: + label - the DMLabel
975: . is - the point IS
976: - value - The point value
978: Level: intermediate
980: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
981: @*/
982: PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value)
983: {
984: PetscInt v, n, p;
985: const PetscInt *points;
986: PetscErrorCode ierr;
991: /* Find label value, add new entry if needed */
992: if (value == label->defaultValue) return(0);
993: DMLabelLookupAddStratum(label, value, &v);
994: /* Set keys */
995: DMLabelMakeInvalid_Private(label, v);
996: ISGetLocalSize(is, &n);
997: ISGetIndices(is, &points);
998: for (p = 0; p < n; ++p) {PetscHSetIAdd(label->ht[v], points[p]);}
999: ISRestoreIndices(is, &points);
1000: return(0);
1001: }
1003: /*@
1004: DMLabelGetNumValues - Get the number of values that the DMLabel takes
1006: Not collective
1008: Input Parameter:
1009: . label - the DMLabel
1011: Output Parameter:
1012: . numValues - the number of values
1014: Level: intermediate
1016: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1017: @*/
1018: PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
1019: {
1023: *numValues = label->numStrata;
1024: return(0);
1025: }
1027: /*@
1028: DMLabelGetValueIS - Get an IS of all values that the DMlabel takes
1030: Not collective
1032: Input Parameter:
1033: . label - the DMLabel
1035: Output Parameter:
1036: . is - the value IS
1038: Level: intermediate
1040: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1041: @*/
1042: PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
1043: {
1049: ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);
1050: return(0);
1051: }
1053: /*@
1054: DMLabelGetValueIndex - Get the index of a given value in the list of values for the DMlabel, or -1 if it is not present
1056: Not collective
1058: Input Parameters:
1059: + label - the DMLabel
1060: - value - the value
1062: Output Parameter:
1063: . index - the index of value in the list of values
1065: Level: intermediate
1067: .seealso: DMLabelGetValueIS(), DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1068: @*/
1069: PetscErrorCode DMLabelGetValueIndex(DMLabel label, PetscInt value, PetscInt *index)
1070: {
1071: PetscInt v;
1076: /* Do not assume they are sorted */
1077: for (v = 0; v < label->numStrata; ++v) if (label->stratumValues[v] == value) break;
1078: if (v >= label->numStrata) *index = -1;
1079: else *index = v;
1080: return(0);
1081: }
1083: /*@
1084: DMLabelHasStratum - Determine whether points exist with the given value
1086: Not collective
1088: Input Parameters:
1089: + label - the DMLabel
1090: - value - the stratum value
1092: Output Parameter:
1093: . exists - Flag saying whether points exist
1095: Level: intermediate
1097: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1098: @*/
1099: PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists)
1100: {
1101: PetscInt v;
1107: DMLabelLookupStratum(label, value, &v);
1108: *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE;
1109: return(0);
1110: }
1112: /*@
1113: DMLabelGetStratumSize - Get the size of a stratum
1115: Not collective
1117: Input Parameters:
1118: + label - the DMLabel
1119: - value - the stratum value
1121: Output Parameter:
1122: . size - The number of points in the stratum
1124: Level: intermediate
1126: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1127: @*/
1128: PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
1129: {
1130: PetscInt v;
1136: *size = 0;
1137: DMLabelLookupStratum(label, value, &v);
1138: if (v < 0) return(0);
1139: DMLabelMakeValid_Private(label, v);
1140: *size = label->stratumSizes[v];
1141: return(0);
1142: }
1144: /*@
1145: DMLabelGetStratumBounds - Get the largest and smallest point of a stratum
1147: Not collective
1149: Input Parameters:
1150: + label - the DMLabel
1151: - value - the stratum value
1153: Output Parameters:
1154: + start - the smallest point in the stratum
1155: - end - the largest point in the stratum
1157: Level: intermediate
1159: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1160: @*/
1161: PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
1162: {
1163: PetscInt v, min, max;
1170: DMLabelLookupStratum(label, value, &v);
1171: if (v < 0) return(0);
1172: DMLabelMakeValid_Private(label, v);
1173: if (label->stratumSizes[v] <= 0) return(0);
1174: ISGetMinMax(label->points[v], &min, &max);
1175: if (start) *start = min;
1176: if (end) *end = max+1;
1177: return(0);
1178: }
1180: /*@
1181: DMLabelGetStratumIS - Get an IS with the stratum points
1183: Not collective
1185: Input Parameters:
1186: + label - the DMLabel
1187: - value - the stratum value
1189: Output Parameter:
1190: . points - The stratum points
1192: Level: intermediate
1194: Notes:
1195: The output IS should be destroyed when no longer needed.
1196: Returns NULL if the stratum is empty.
1198: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1199: @*/
1200: PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
1201: {
1202: PetscInt v;
1208: *points = NULL;
1209: DMLabelLookupStratum(label, value, &v);
1210: if (v < 0) return(0);
1211: DMLabelMakeValid_Private(label, v);
1212: PetscObjectReference((PetscObject) label->points[v]);
1213: *points = label->points[v];
1214: return(0);
1215: }
1217: /*@
1218: DMLabelSetStratumIS - Set the stratum points using an IS
1220: Not collective
1222: Input Parameters:
1223: + label - the DMLabel
1224: . value - the stratum value
1225: - points - The stratum points
1227: Level: intermediate
1229: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1230: @*/
1231: PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is)
1232: {
1233: PetscInt v;
1239: DMLabelLookupAddStratum(label, value, &v);
1240: if (is == label->points[v]) return(0);
1241: DMLabelClearStratum(label, value);
1242: ISGetLocalSize(is, &(label->stratumSizes[v]));
1243: PetscObjectReference((PetscObject)is);
1244: ISDestroy(&(label->points[v]));
1245: label->points[v] = is;
1246: label->validIS[v] = PETSC_TRUE;
1247: PetscObjectStateIncrease((PetscObject) label);
1248: if (label->bt) {
1249: const PetscInt *points;
1250: PetscInt p;
1252: ISGetIndices(is,&points);
1253: for (p = 0; p < label->stratumSizes[v]; ++p) {
1254: const PetscInt point = points[p];
1256: 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);
1257: PetscBTSet(label->bt, point - label->pStart);
1258: }
1259: }
1260: return(0);
1261: }
1263: /*@
1264: DMLabelClearStratum - Remove a stratum
1266: Not collective
1268: Input Parameters:
1269: + label - the DMLabel
1270: - value - the stratum value
1272: Level: intermediate
1274: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1275: @*/
1276: PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
1277: {
1278: PetscInt v;
1283: DMLabelLookupStratum(label, value, &v);
1284: if (v < 0) return(0);
1285: if (label->validIS[v]) {
1286: if (label->bt) {
1287: PetscInt i;
1288: const PetscInt *points;
1290: ISGetIndices(label->points[v], &points);
1291: for (i = 0; i < label->stratumSizes[v]; ++i) {
1292: const PetscInt point = points[i];
1294: 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);
1295: PetscBTClear(label->bt, point - label->pStart);
1296: }
1297: ISRestoreIndices(label->points[v], &points);
1298: }
1299: label->stratumSizes[v] = 0;
1300: ISDestroy(&label->points[v]);
1301: ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &label->points[v]);
1302: PetscObjectSetName((PetscObject) label->points[v], "indices");
1303: PetscObjectStateIncrease((PetscObject) label);
1304: } else {
1305: PetscHSetIClear(label->ht[v]);
1306: }
1307: return(0);
1308: }
1310: /*@
1311: DMLabelSetStratumBounds - Efficiently give a contiguous set of points a given label value
1313: Not collective
1315: Input Parameters:
1316: + label - The DMLabel
1317: . value - The label value for all points
1318: . pStart - The first point
1319: - pEnd - A point beyond all marked points
1321: Note: The marks points are [pStart, pEnd), and only the bounds are stored.
1323: Level: intermediate
1325: .seealso: DMLabelCreate(), DMLabelSetStratumIS(), DMLabelGetStratumIS()
1326: @*/
1327: PetscErrorCode DMLabelSetStratumBounds(DMLabel label, PetscInt value, PetscInt pStart, PetscInt pEnd)
1328: {
1329: IS pIS;
1333: ISCreateStride(PETSC_COMM_SELF, pEnd - pStart, pStart, 1, &pIS);
1334: DMLabelSetStratumIS(label, value, pIS);
1335: ISDestroy(&pIS);
1336: return(0);
1337: }
1339: /*@
1340: DMLabelGetStratumPointIndex - Get the index of a point in a given stratum
1342: Not collective
1344: Input Parameters:
1345: + label - The DMLabel
1346: . value - The label value
1347: - p - A point with this value
1349: Output Parameter:
1350: . index - The index of this point in the stratum, or -1 if the point is not in the stratum or the stratum does not exist
1352: Level: intermediate
1354: .seealso: DMLabelGetValueIndex(), DMLabelGetStratumIS(), DMLabelCreate()
1355: @*/
1356: PetscErrorCode DMLabelGetStratumPointIndex(DMLabel label, PetscInt value, PetscInt p, PetscInt *index)
1357: {
1358: const PetscInt *indices;
1359: PetscInt v;
1360: PetscErrorCode ierr;
1365: *index = -1;
1366: DMLabelLookupStratum(label, value, &v);
1367: if (v < 0) return(0);
1368: DMLabelMakeValid_Private(label, v);
1369: ISGetIndices(label->points[v], &indices);
1370: PetscFindInt(p, label->stratumSizes[v], indices, index);
1371: ISRestoreIndices(label->points[v], &indices);
1372: return(0);
1373: }
1375: /*@
1376: DMLabelFilter - Remove all points outside of [start, end)
1378: Not collective
1380: Input Parameters:
1381: + label - the DMLabel
1382: . start - the first point kept
1383: - end - one more than the last point kept
1385: Level: intermediate
1387: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1388: @*/
1389: PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
1390: {
1391: PetscInt v;
1396: DMLabelDestroyIndex(label);
1397: DMLabelMakeAllValid_Private(label);
1398: for (v = 0; v < label->numStrata; ++v) {
1399: ISGeneralFilter(label->points[v], start, end);
1400: }
1401: DMLabelCreateIndex(label, start, end);
1402: return(0);
1403: }
1405: /*@
1406: DMLabelPermute - Create a new label with permuted points
1408: Not collective
1410: Input Parameters:
1411: + label - the DMLabel
1412: - permutation - the point permutation
1414: Output Parameter:
1415: . labelnew - the new label containing the permuted points
1417: Level: intermediate
1419: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1420: @*/
1421: PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
1422: {
1423: const PetscInt *perm;
1424: PetscInt numValues, numPoints, v, q;
1425: PetscErrorCode ierr;
1430: DMLabelMakeAllValid_Private(label);
1431: DMLabelDuplicate(label, labelNew);
1432: DMLabelGetNumValues(*labelNew, &numValues);
1433: ISGetLocalSize(permutation, &numPoints);
1434: ISGetIndices(permutation, &perm);
1435: for (v = 0; v < numValues; ++v) {
1436: const PetscInt size = (*labelNew)->stratumSizes[v];
1437: const PetscInt *points;
1438: PetscInt *pointsNew;
1440: ISGetIndices((*labelNew)->points[v],&points);
1441: PetscMalloc1(size,&pointsNew);
1442: for (q = 0; q < size; ++q) {
1443: const PetscInt point = points[q];
1445: 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);
1446: pointsNew[q] = perm[point];
1447: }
1448: ISRestoreIndices((*labelNew)->points[v],&points);
1449: PetscSortInt(size, pointsNew);
1450: ISDestroy(&((*labelNew)->points[v]));
1451: if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) {
1452: ISCreateStride(PETSC_COMM_SELF,size,pointsNew[0],1,&((*labelNew)->points[v]));
1453: PetscFree(pointsNew);
1454: } else {
1455: ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));
1456: }
1457: PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");
1458: }
1459: ISRestoreIndices(permutation, &perm);
1460: if (label->bt) {
1461: PetscBTDestroy(&label->bt);
1462: DMLabelCreateIndex(label, label->pStart, label->pEnd);
1463: }
1464: return(0);
1465: }
1467: PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata)
1468: {
1469: MPI_Comm comm;
1470: PetscInt s, l, nroots, nleaves, dof, offset, size;
1471: PetscInt *remoteOffsets, *rootStrata, *rootIdx;
1472: PetscSection rootSection;
1473: PetscSF labelSF;
1477: if (label) {DMLabelMakeAllValid_Private(label);}
1478: PetscObjectGetComm((PetscObject)sf, &comm);
1479: /* Build a section of stratum values per point, generate the according SF
1480: and distribute point-wise stratum values to leaves. */
1481: PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);
1482: PetscSectionCreate(comm, &rootSection);
1483: PetscSectionSetChart(rootSection, 0, nroots);
1484: if (label) {
1485: for (s = 0; s < label->numStrata; ++s) {
1486: const PetscInt *points;
1488: ISGetIndices(label->points[s], &points);
1489: for (l = 0; l < label->stratumSizes[s]; l++) {
1490: PetscSectionGetDof(rootSection, points[l], &dof);
1491: PetscSectionSetDof(rootSection, points[l], dof+1);
1492: }
1493: ISRestoreIndices(label->points[s], &points);
1494: }
1495: }
1496: PetscSectionSetUp(rootSection);
1497: /* Create a point-wise array of stratum values */
1498: PetscSectionGetStorageSize(rootSection, &size);
1499: PetscMalloc1(size, &rootStrata);
1500: PetscCalloc1(nroots, &rootIdx);
1501: if (label) {
1502: for (s = 0; s < label->numStrata; ++s) {
1503: const PetscInt *points;
1505: ISGetIndices(label->points[s], &points);
1506: for (l = 0; l < label->stratumSizes[s]; l++) {
1507: const PetscInt p = points[l];
1508: PetscSectionGetOffset(rootSection, p, &offset);
1509: rootStrata[offset+rootIdx[p]++] = label->stratumValues[s];
1510: }
1511: ISRestoreIndices(label->points[s], &points);
1512: }
1513: }
1514: /* Build SF that maps label points to remote processes */
1515: PetscSectionCreate(comm, leafSection);
1516: PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);
1517: PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);
1518: PetscFree(remoteOffsets);
1519: /* Send the strata for each point over the derived SF */
1520: PetscSectionGetStorageSize(*leafSection, &size);
1521: PetscMalloc1(size, leafStrata);
1522: PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata,MPI_REPLACE);
1523: PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata,MPI_REPLACE);
1524: /* Clean up */
1525: PetscFree(rootStrata);
1526: PetscFree(rootIdx);
1527: PetscSectionDestroy(&rootSection);
1528: PetscSFDestroy(&labelSF);
1529: return(0);
1530: }
1532: /*@
1533: DMLabelDistribute - Create a new label pushed forward over the PetscSF
1535: Collective on sf
1537: Input Parameters:
1538: + label - the DMLabel
1539: - sf - the map from old to new distribution
1541: Output Parameter:
1542: . labelnew - the new redistributed label
1544: Level: intermediate
1546: .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1547: @*/
1548: PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew)
1549: {
1550: MPI_Comm comm;
1551: PetscSection leafSection;
1552: PetscInt p, pStart, pEnd, s, size, dof, offset, stratum;
1553: PetscInt *leafStrata, *strataIdx;
1554: PetscInt **points;
1555: const char *lname = NULL;
1556: char *name;
1557: PetscInt nameSize;
1558: PetscHSetI stratumHash;
1559: size_t len = 0;
1560: PetscMPIInt rank;
1565: if (label) {
1567: DMLabelMakeAllValid_Private(label);
1568: }
1569: PetscObjectGetComm((PetscObject)sf, &comm);
1570: MPI_Comm_rank(comm, &rank);
1571: /* Bcast name */
1572: if (rank == 0) {
1573: PetscObjectGetName((PetscObject) label, &lname);
1574: PetscStrlen(lname, &len);
1575: }
1576: nameSize = len;
1577: MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);
1578: PetscMalloc1(nameSize+1, &name);
1579: if (rank == 0) {PetscArraycpy(name, lname, nameSize+1);}
1580: MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);
1581: DMLabelCreate(PETSC_COMM_SELF, name, labelNew);
1582: PetscFree(name);
1583: /* Bcast defaultValue */
1584: if (rank == 0) (*labelNew)->defaultValue = label->defaultValue;
1585: MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);
1586: /* Distribute stratum values over the SF and get the point mapping on the receiver */
1587: DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);
1588: /* Determine received stratum values and initialise new label*/
1589: PetscHSetICreate(&stratumHash);
1590: PetscSectionGetStorageSize(leafSection, &size);
1591: for (p = 0; p < size; ++p) {PetscHSetIAdd(stratumHash, leafStrata[p]);}
1592: PetscHSetIGetSize(stratumHash, &(*labelNew)->numStrata);
1593: PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);
1594: for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE;
1595: PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);
1596: /* Turn leafStrata into indices rather than stratum values */
1597: offset = 0;
1598: PetscHSetIGetElems(stratumHash, &offset, (*labelNew)->stratumValues);
1599: PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);
1600: for (s = 0; s < (*labelNew)->numStrata; ++s) {
1601: PetscHMapISet((*labelNew)->hmap, (*labelNew)->stratumValues[s], s);
1602: }
1603: for (p = 0; p < size; ++p) {
1604: for (s = 0; s < (*labelNew)->numStrata; ++s) {
1605: if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;}
1606: }
1607: }
1608: /* Rebuild the point strata on the receiver */
1609: PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);
1610: PetscSectionGetChart(leafSection, &pStart, &pEnd);
1611: for (p=pStart; p<pEnd; p++) {
1612: PetscSectionGetDof(leafSection, p, &dof);
1613: PetscSectionGetOffset(leafSection, p, &offset);
1614: for (s=0; s<dof; s++) {
1615: (*labelNew)->stratumSizes[leafStrata[offset+s]]++;
1616: }
1617: }
1618: PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);
1619: PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);
1620: PetscMalloc1((*labelNew)->numStrata,&points);
1621: for (s = 0; s < (*labelNew)->numStrata; ++s) {
1622: PetscHSetICreate(&(*labelNew)->ht[s]);
1623: PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));
1624: }
1625: /* Insert points into new strata */
1626: PetscCalloc1((*labelNew)->numStrata, &strataIdx);
1627: PetscSectionGetChart(leafSection, &pStart, &pEnd);
1628: for (p=pStart; p<pEnd; p++) {
1629: PetscSectionGetDof(leafSection, p, &dof);
1630: PetscSectionGetOffset(leafSection, p, &offset);
1631: for (s=0; s<dof; s++) {
1632: stratum = leafStrata[offset+s];
1633: points[stratum][strataIdx[stratum]++] = p;
1634: }
1635: }
1636: for (s = 0; s < (*labelNew)->numStrata; s++) {
1637: ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));
1638: PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");
1639: }
1640: PetscFree(points);
1641: PetscHSetIDestroy(&stratumHash);
1642: PetscFree(leafStrata);
1643: PetscFree(strataIdx);
1644: PetscSectionDestroy(&leafSection);
1645: return(0);
1646: }
1648: /*@
1649: DMLabelGather - Gather all label values from leafs into roots
1651: Collective on sf
1653: Input Parameters:
1654: + label - the DMLabel
1655: - sf - the Star Forest point communication map
1657: Output Parameters:
1658: . labelNew - the new DMLabel with localised leaf values
1660: Level: developer
1662: Note: This is the inverse operation to DMLabelDistribute.
1664: .seealso: DMLabelDistribute()
1665: @*/
1666: PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew)
1667: {
1668: MPI_Comm comm;
1669: PetscSection rootSection;
1670: PetscSF sfLabel;
1671: PetscSFNode *rootPoints, *leafPoints;
1672: PetscInt p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset;
1673: const PetscInt *rootDegree, *ilocal;
1674: PetscInt *rootStrata;
1675: const char *lname;
1676: char *name;
1677: PetscInt nameSize;
1678: size_t len = 0;
1679: PetscMPIInt rank, size;
1685: PetscObjectGetComm((PetscObject)sf, &comm);
1686: MPI_Comm_rank(comm, &rank);
1687: MPI_Comm_size(comm, &size);
1688: /* Bcast name */
1689: if (rank == 0) {
1690: PetscObjectGetName((PetscObject) label, &lname);
1691: PetscStrlen(lname, &len);
1692: }
1693: nameSize = len;
1694: MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);
1695: PetscMalloc1(nameSize+1, &name);
1696: if (rank == 0) {PetscArraycpy(name, lname, nameSize+1);}
1697: MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);
1698: DMLabelCreate(PETSC_COMM_SELF, name, labelNew);
1699: PetscFree(name);
1700: /* Gather rank/index pairs of leaves into local roots to build
1701: an inverse, multi-rooted SF. Note that this ignores local leaf
1702: indexing due to the use of the multiSF in PetscSFGather. */
1703: PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);
1704: PetscMalloc1(nroots, &leafPoints);
1705: for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1;
1706: for (p = 0; p < nleaves; p++) {
1707: PetscInt ilp = ilocal ? ilocal[p] : p;
1709: leafPoints[ilp].index = ilp;
1710: leafPoints[ilp].rank = rank;
1711: }
1712: PetscSFComputeDegreeBegin(sf, &rootDegree);
1713: PetscSFComputeDegreeEnd(sf, &rootDegree);
1714: for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p];
1715: PetscMalloc1(nmultiroots, &rootPoints);
1716: PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);
1717: PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);
1718: PetscSFCreate(comm,& sfLabel);
1719: PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);
1720: /* Migrate label over inverted SF to pull stratum values at leaves into roots. */
1721: DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);
1722: /* Rebuild the point strata on the receiver */
1723: for (p = 0, idx = 0; p < nroots; p++) {
1724: for (d = 0; d < rootDegree[p]; d++) {
1725: PetscSectionGetDof(rootSection, idx+d, &dof);
1726: PetscSectionGetOffset(rootSection, idx+d, &offset);
1727: for (s = 0; s < dof; s++) {DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);}
1728: }
1729: idx += rootDegree[p];
1730: }
1731: PetscFree(leafPoints);
1732: PetscFree(rootStrata);
1733: PetscSectionDestroy(&rootSection);
1734: PetscSFDestroy(&sfLabel);
1735: return(0);
1736: }
1738: /*@
1739: DMLabelConvertToSection - Make a PetscSection/IS pair that encodes the label
1741: Not collective
1743: Input Parameter:
1744: . label - the DMLabel
1746: Output Parameters:
1747: + section - the section giving offsets for each stratum
1748: - is - An IS containing all the label points
1750: Level: developer
1752: .seealso: DMLabelDistribute()
1753: @*/
1754: PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is)
1755: {
1756: IS vIS;
1757: const PetscInt *values;
1758: PetscInt *points;
1759: PetscInt nV, vS = 0, vE = 0, v, N;
1760: PetscErrorCode ierr;
1764: DMLabelGetNumValues(label, &nV);
1765: DMLabelGetValueIS(label, &vIS);
1766: ISGetIndices(vIS, &values);
1767: if (nV) {vS = values[0]; vE = values[0]+1;}
1768: for (v = 1; v < nV; ++v) {
1769: vS = PetscMin(vS, values[v]);
1770: vE = PetscMax(vE, values[v]+1);
1771: }
1772: PetscSectionCreate(PETSC_COMM_SELF, section);
1773: PetscSectionSetChart(*section, vS, vE);
1774: for (v = 0; v < nV; ++v) {
1775: PetscInt n;
1777: DMLabelGetStratumSize(label, values[v], &n);
1778: PetscSectionSetDof(*section, values[v], n);
1779: }
1780: PetscSectionSetUp(*section);
1781: PetscSectionGetStorageSize(*section, &N);
1782: PetscMalloc1(N, &points);
1783: for (v = 0; v < nV; ++v) {
1784: IS is;
1785: const PetscInt *spoints;
1786: PetscInt dof, off, p;
1788: PetscSectionGetDof(*section, values[v], &dof);
1789: PetscSectionGetOffset(*section, values[v], &off);
1790: DMLabelGetStratumIS(label, values[v], &is);
1791: ISGetIndices(is, &spoints);
1792: for (p = 0; p < dof; ++p) points[off+p] = spoints[p];
1793: ISRestoreIndices(is, &spoints);
1794: ISDestroy(&is);
1795: }
1796: ISRestoreIndices(vIS, &values);
1797: ISDestroy(&vIS);
1798: ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);
1799: return(0);
1800: }
1802: /*@
1803: PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using
1804: the local section and an SF describing the section point overlap.
1806: Collective on sf
1808: Input Parameters:
1809: + s - The PetscSection for the local field layout
1810: . sf - The SF describing parallel layout of the section points
1811: . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1812: . label - The label specifying the points
1813: - labelValue - The label stratum specifying the points
1815: Output Parameter:
1816: . gsection - The PetscSection for the global field layout
1818: Note: This gives negative sizes and offsets to points not owned by this process
1820: Level: developer
1822: .seealso: PetscSectionCreate()
1823: @*/
1824: PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection)
1825: {
1826: PetscInt *neg = NULL, *tmpOff = NULL;
1827: PetscInt pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots;
1834: PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);
1835: PetscSectionGetChart(s, &pStart, &pEnd);
1836: PetscSectionSetChart(*gsection, pStart, pEnd);
1837: PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);
1838: if (nroots >= 0) {
1839: if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart);
1840: PetscCalloc1(nroots, &neg);
1841: if (nroots > pEnd-pStart) {
1842: PetscCalloc1(nroots, &tmpOff);
1843: } else {
1844: tmpOff = &(*gsection)->atlasDof[-pStart];
1845: }
1846: }
1847: /* Mark ghost points with negative dof */
1848: for (p = pStart; p < pEnd; ++p) {
1849: PetscInt value;
1851: DMLabelGetValue(label, p, &value);
1852: if (value != labelValue) continue;
1853: PetscSectionGetDof(s, p, &dof);
1854: PetscSectionSetDof(*gsection, p, dof);
1855: PetscSectionGetConstraintDof(s, p, &cdof);
1856: if (!includeConstraints && cdof > 0) {PetscSectionSetConstraintDof(*gsection, p, cdof);}
1857: if (neg) neg[p] = -(dof+1);
1858: }
1859: PetscSectionSetUpBC(*gsection);
1860: if (nroots >= 0) {
1861: PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff,MPI_REPLACE);
1862: PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff,MPI_REPLACE);
1863: if (nroots > pEnd-pStart) {
1864: for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
1865: }
1866: }
1867: /* Calculate new sizes, get proccess offset, and calculate point offsets */
1868: for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1869: cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
1870: (*gsection)->atlasOff[p] = off;
1871: off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
1872: }
1873: MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));
1874: globalOff -= off;
1875: for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1876: (*gsection)->atlasOff[p] += globalOff;
1877: if (neg) neg[p] = -((*gsection)->atlasOff[p]+1);
1878: }
1879: /* Put in negative offsets for ghost points */
1880: if (nroots >= 0) {
1881: PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff,MPI_REPLACE);
1882: PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff,MPI_REPLACE);
1883: if (nroots > pEnd-pStart) {
1884: for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1885: }
1886: }
1887: if (nroots >= 0 && nroots > pEnd-pStart) {PetscFree(tmpOff);}
1888: PetscFree(neg);
1889: return(0);
1890: }
1892: typedef struct _n_PetscSectionSym_Label
1893: {
1894: DMLabel label;
1895: PetscCopyMode *modes;
1896: PetscInt *sizes;
1897: const PetscInt ***perms;
1898: const PetscScalar ***rots;
1899: PetscInt (*minMaxOrients)[2];
1900: PetscInt numStrata; /* numStrata is only increasing, functions as a state */
1901: } PetscSectionSym_Label;
1903: static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym)
1904: {
1905: PetscInt i, j;
1906: PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
1907: PetscErrorCode ierr;
1910: for (i = 0; i <= sl->numStrata; i++) {
1911: if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) {
1912: for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
1913: if (sl->perms[i]) {PetscFree(sl->perms[i][j]);}
1914: if (sl->rots[i]) {PetscFree(sl->rots[i][j]);}
1915: }
1916: if (sl->perms[i]) {
1917: const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]];
1919: PetscFree(perms);
1920: }
1921: if (sl->rots[i]) {
1922: const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]];
1924: PetscFree(rots);
1925: }
1926: }
1927: }
1928: PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);
1929: DMLabelDestroy(&sl->label);
1930: sl->numStrata = 0;
1931: return(0);
1932: }
1934: static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym)
1935: {
1939: PetscSectionSymLabelReset(sym);
1940: PetscFree(sym->data);
1941: return(0);
1942: }
1944: static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer)
1945: {
1946: PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
1947: PetscBool isAscii;
1948: DMLabel label = sl->label;
1949: const char *name;
1950: PetscErrorCode ierr;
1953: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);
1954: if (isAscii) {
1955: PetscInt i, j, k;
1956: PetscViewerFormat format;
1958: PetscViewerGetFormat(viewer,&format);
1959: if (label) {
1960: PetscViewerGetFormat(viewer,&format);
1961: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1962: PetscViewerASCIIPushTab(viewer);
1963: DMLabelView(label, viewer);
1964: PetscViewerASCIIPopTab(viewer);
1965: } else {
1966: PetscObjectGetName((PetscObject) sl->label, &name);
1967: PetscViewerASCIIPrintf(viewer," Label '%s'\n",name);
1968: }
1969: } else {
1970: PetscViewerASCIIPrintf(viewer, "No label given\n");
1971: }
1972: PetscViewerASCIIPushTab(viewer);
1973: for (i = 0; i <= sl->numStrata; i++) {
1974: PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue;
1976: if (!(sl->perms[i] || sl->rots[i])) {
1977: PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);
1978: } else {
1979: PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);
1980: PetscViewerASCIIPushTab(viewer);
1981: PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);
1982: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1983: PetscViewerASCIIPushTab(viewer);
1984: for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
1985: if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) {
1986: PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);
1987: } else {
1988: PetscInt tab;
1990: PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);
1991: PetscViewerASCIIPushTab(viewer);
1992: PetscViewerASCIIGetTab(viewer,&tab);
1993: if (sl->perms[i] && sl->perms[i][j]) {
1994: PetscViewerASCIIPrintf(viewer,"Permutation:");
1995: PetscViewerASCIISetTab(viewer,0);
1996: for (k = 0; k < sl->sizes[i]; k++) {PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);}
1997: PetscViewerASCIIPrintf(viewer,"\n");
1998: PetscViewerASCIISetTab(viewer,tab);
1999: }
2000: if (sl->rots[i] && sl->rots[i][j]) {
2001: PetscViewerASCIIPrintf(viewer,"Rotations: ");
2002: PetscViewerASCIISetTab(viewer,0);
2003: #if defined(PETSC_USE_COMPLEX)
2004: 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]));}
2005: #else
2006: for (k = 0; k < sl->sizes[i]; k++) {PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);}
2007: #endif
2008: PetscViewerASCIIPrintf(viewer,"\n");
2009: PetscViewerASCIISetTab(viewer,tab);
2010: }
2011: PetscViewerASCIIPopTab(viewer);
2012: }
2013: }
2014: PetscViewerASCIIPopTab(viewer);
2015: }
2016: PetscViewerASCIIPopTab(viewer);
2017: }
2018: }
2019: PetscViewerASCIIPopTab(viewer);
2020: }
2021: return(0);
2022: }
2024: /*@
2025: PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries
2027: Logically collective on sym
2029: Input parameters:
2030: + sym - the section symmetries
2031: - label - the DMLabel describing the types of points
2033: Level: developer:
2035: .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms()
2036: @*/
2037: PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label)
2038: {
2039: PetscSectionSym_Label *sl;
2040: PetscErrorCode ierr;
2044: sl = (PetscSectionSym_Label *) sym->data;
2045: if (sl->label && sl->label != label) {PetscSectionSymLabelReset(sym);}
2046: if (label) {
2047: sl->label = label;
2048: PetscObjectReference((PetscObject) label);
2049: DMLabelGetNumValues(label,&sl->numStrata);
2050: 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);
2051: PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));
2052: PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));
2053: PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));
2054: PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));
2055: PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));
2056: }
2057: return(0);
2058: }
2060: /*@C
2061: PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum
2063: Logically collective on sym
2065: InputParameters:
2066: + sym - the section symmetries
2067: . stratum - the stratum value in the label that we are assigning symmetries for
2068: . size - the number of dofs for points in the stratum of the label
2069: . minOrient - the smallest orientation for a point in this stratum
2070: . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient))
2071: . mode - how sym should copy the perms and rots arrays
2072: . perms - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation. A NULL permutation is the identity
2073: - 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
2075: Level: developer
2077: .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel()
2078: @*/
2079: PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots)
2080: {
2081: PetscSectionSym_Label *sl;
2082: const char *name;
2083: PetscInt i, j, k;
2084: PetscErrorCode ierr;
2088: sl = (PetscSectionSym_Label *) sym->data;
2089: if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet");
2090: for (i = 0; i <= sl->numStrata; i++) {
2091: PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue;
2093: if (stratum == value) break;
2094: }
2095: PetscObjectGetName((PetscObject) sl->label, &name);
2096: if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,name);
2097: sl->sizes[i] = size;
2098: sl->modes[i] = mode;
2099: sl->minMaxOrients[i][0] = minOrient;
2100: sl->minMaxOrients[i][1] = maxOrient;
2101: if (mode == PETSC_COPY_VALUES) {
2102: if (perms) {
2103: PetscInt **ownPerms;
2105: PetscCalloc1(maxOrient - minOrient,&ownPerms);
2106: for (j = 0; j < maxOrient-minOrient; j++) {
2107: if (perms[j]) {
2108: PetscMalloc1(size,&ownPerms[j]);
2109: for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];}
2110: }
2111: }
2112: sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient];
2113: }
2114: if (rots) {
2115: PetscScalar **ownRots;
2117: PetscCalloc1(maxOrient - minOrient,&ownRots);
2118: for (j = 0; j < maxOrient-minOrient; j++) {
2119: if (rots[j]) {
2120: PetscMalloc1(size,&ownRots[j]);
2121: for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];}
2122: }
2123: }
2124: sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient];
2125: }
2126: } else {
2127: sl->perms[i] = perms ? &perms[-minOrient] : NULL;
2128: sl->rots[i] = rots ? &rots[-minOrient] : NULL;
2129: }
2130: return(0);
2131: }
2133: static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots)
2134: {
2135: PetscInt i, j, numStrata;
2136: PetscSectionSym_Label *sl;
2137: DMLabel label;
2138: PetscErrorCode ierr;
2141: sl = (PetscSectionSym_Label *) sym->data;
2142: numStrata = sl->numStrata;
2143: label = sl->label;
2144: for (i = 0; i < numPoints; i++) {
2145: PetscInt point = points[2*i];
2146: PetscInt ornt = points[2*i+1];
2148: for (j = 0; j < numStrata; j++) {
2149: if (label->validIS[j]) {
2150: PetscInt k;
2152: ISLocate(label->points[j],point,&k);
2153: if (k >= 0) break;
2154: } else {
2155: PetscBool has;
2157: PetscHSetIHas(label->ht[j], point, &has);
2158: if (has) break;
2159: }
2160: }
2161: 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);
2162: if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;}
2163: if (rots) {rots[i] = sl->rots[j] ? sl->rots[j][ornt] : NULL;}
2164: }
2165: return(0);
2166: }
2168: PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym)
2169: {
2170: PetscSectionSym_Label *sl;
2171: PetscErrorCode ierr;
2174: PetscNewLog(sym,&sl);
2175: sym->ops->getpoints = PetscSectionSymGetPoints_Label;
2176: sym->ops->view = PetscSectionSymView_Label;
2177: sym->ops->destroy = PetscSectionSymDestroy_Label;
2178: sym->data = (void *) sl;
2179: return(0);
2180: }
2182: /*@
2183: PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label
2185: Collective
2187: Input Parameters:
2188: + comm - the MPI communicator for the new symmetry
2189: - label - the label defining the strata
2191: Output Parameters:
2192: . sym - the section symmetries
2194: Level: developer
2196: .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()
2197: @*/
2198: PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym)
2199: {
2200: PetscErrorCode ierr;
2203: DMInitializePackage();
2204: PetscSectionSymCreate(comm,sym);
2205: PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);
2206: PetscSectionSymLabelSetLabel(*sym,label);
2207: return(0);
2208: }