19f6c5813SMatthew G. Knepley #include <petsc/private/dmlabelimpl.h> /*I "petscdmlabelephemeral.h" I*/
29f6c5813SMatthew G. Knepley #include <petscdmlabelephemeral.h>
39f6c5813SMatthew G. Knepley
49f6c5813SMatthew G. Knepley /*
59f6c5813SMatthew G. Knepley An emphemeral label is read-only
69f6c5813SMatthew G. Knepley
79f6c5813SMatthew G. Knepley This initial implementation makes the assumption that the produced points have unique parents. If this is
89f6c5813SMatthew G. Knepley not satisfied, hash tables can be introduced to ensure produced points are unique.
99f6c5813SMatthew G. Knepley */
108dcf10e8SMatthew G. Knepley static PetscErrorCode DMLabelInitialize_Ephemeral(DMLabel);
119f6c5813SMatthew G. Knepley
DMLabelEphemeralComputeStratumSize_Private(DMLabel label,PetscInt value)129f6c5813SMatthew G. Knepley static PetscErrorCode DMLabelEphemeralComputeStratumSize_Private(DMLabel label, PetscInt value)
139f6c5813SMatthew G. Knepley {
149f6c5813SMatthew G. Knepley DMPlexTransform tr;
159f6c5813SMatthew G. Knepley DM dm;
169f6c5813SMatthew G. Knepley DMLabel olabel;
179f6c5813SMatthew G. Knepley IS opointIS;
189f6c5813SMatthew G. Knepley const PetscInt *opoints;
199f6c5813SMatthew G. Knepley PetscInt Np = 0, Nop, op, v;
209f6c5813SMatthew G. Knepley
219f6c5813SMatthew G. Knepley PetscFunctionBegin;
229f6c5813SMatthew G. Knepley PetscCall(DMLabelEphemeralGetTransform(label, &tr));
239f6c5813SMatthew G. Knepley PetscCall(DMLabelEphemeralGetLabel(label, &olabel));
249f6c5813SMatthew G. Knepley PetscCall(DMPlexTransformGetDM(tr, &dm));
259f6c5813SMatthew G. Knepley
269f6c5813SMatthew G. Knepley PetscCall(DMLabelLookupStratum(olabel, value, &v));
279f6c5813SMatthew G. Knepley PetscCall(DMLabelGetStratumIS(olabel, value, &opointIS));
289f6c5813SMatthew G. Knepley PetscCall(ISGetLocalSize(opointIS, &Nop));
299f6c5813SMatthew G. Knepley PetscCall(ISGetIndices(opointIS, &opoints));
309f6c5813SMatthew G. Knepley for (op = 0; op < Nop; ++op) {
319f6c5813SMatthew G. Knepley const PetscInt point = opoints[op];
329f6c5813SMatthew G. Knepley DMPolytopeType ct;
339f6c5813SMatthew G. Knepley DMPolytopeType *rct;
349f6c5813SMatthew G. Knepley PetscInt *rsize, *rcone, *rornt;
359f6c5813SMatthew G. Knepley PetscInt Nct;
369f6c5813SMatthew G. Knepley
379f6c5813SMatthew G. Knepley PetscCall(DMPlexGetCellType(dm, point, &ct));
389f6c5813SMatthew G. Knepley PetscCall(DMPlexTransformCellTransform(tr, ct, point, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
399f6c5813SMatthew G. Knepley for (PetscInt n = 0; n < Nct; ++n) Np += rsize[n];
409f6c5813SMatthew G. Knepley }
419f6c5813SMatthew G. Knepley PetscCall(ISRestoreIndices(opointIS, &opoints));
429f6c5813SMatthew G. Knepley PetscCall(ISDestroy(&opointIS));
439f6c5813SMatthew G. Knepley label->stratumSizes[v] = Np;
443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
459f6c5813SMatthew G. Knepley }
469f6c5813SMatthew G. Knepley
DMLabelGetStratumIS_Ephemeral(DMLabel label,PetscInt v,IS * stratum)479f6c5813SMatthew G. Knepley static PetscErrorCode DMLabelGetStratumIS_Ephemeral(DMLabel label, PetscInt v, IS *stratum)
489f6c5813SMatthew G. Knepley {
499f6c5813SMatthew G. Knepley DMPlexTransform tr;
509f6c5813SMatthew G. Knepley DM dm;
519f6c5813SMatthew G. Knepley DMLabel olabel;
529f6c5813SMatthew G. Knepley IS opointIS;
539f6c5813SMatthew G. Knepley const PetscInt *opoints;
549f6c5813SMatthew G. Knepley PetscInt *points;
559f6c5813SMatthew G. Knepley PetscInt Np, p, Nop, op;
569f6c5813SMatthew G. Knepley
579f6c5813SMatthew G. Knepley PetscFunctionBegin;
589f6c5813SMatthew G. Knepley PetscCall(DMLabelEphemeralGetTransform(label, &tr));
599f6c5813SMatthew G. Knepley PetscCall(DMLabelEphemeralGetLabel(label, &olabel));
609f6c5813SMatthew G. Knepley PetscCall(DMPlexTransformGetDM(tr, &dm));
619f6c5813SMatthew G. Knepley
629f6c5813SMatthew G. Knepley PetscCall(DMLabelGetStratumSize_Private(label, v, &Np));
639f6c5813SMatthew G. Knepley PetscCall(PetscMalloc1(Np, &points));
649f6c5813SMatthew G. Knepley PetscUseTypeMethod(olabel, getstratumis, v, &opointIS);
659f6c5813SMatthew G. Knepley PetscCall(ISGetLocalSize(opointIS, &Nop));
669f6c5813SMatthew G. Knepley PetscCall(ISGetIndices(opointIS, &opoints));
679f6c5813SMatthew G. Knepley for (op = 0, p = 0; op < Nop; ++op) {
689f6c5813SMatthew G. Knepley const PetscInt point = opoints[op];
699f6c5813SMatthew G. Knepley DMPolytopeType ct;
709f6c5813SMatthew G. Knepley DMPolytopeType *rct;
719f6c5813SMatthew G. Knepley PetscInt *rsize, *rcone, *rornt;
729f6c5813SMatthew G. Knepley PetscInt Nct, n, r, pNew = 0;
739f6c5813SMatthew G. Knepley
749f6c5813SMatthew G. Knepley PetscCall(DMPlexGetCellType(dm, point, &ct));
759f6c5813SMatthew G. Knepley PetscCall(DMPlexTransformCellTransform(tr, ct, point, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
769f6c5813SMatthew G. Knepley for (n = 0; n < Nct; ++n) {
779f6c5813SMatthew G. Knepley for (r = 0; r < rsize[n]; ++r) {
789f6c5813SMatthew G. Knepley PetscCall(DMPlexTransformGetTargetPoint(tr, ct, rct[n], point, r, &pNew));
799f6c5813SMatthew G. Knepley points[p++] = pNew;
809f6c5813SMatthew G. Knepley }
819f6c5813SMatthew G. Knepley }
829f6c5813SMatthew G. Knepley }
839f6c5813SMatthew G. Knepley PetscCheck(p == Np, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of stratum points %" PetscInt_FMT " != %" PetscInt_FMT " precomputed size", p, Np);
849f6c5813SMatthew G. Knepley PetscCall(ISRestoreIndices(opointIS, &opoints));
859f6c5813SMatthew G. Knepley PetscCall(ISDestroy(&opointIS));
869f6c5813SMatthew G. Knepley PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Np, points, PETSC_OWN_POINTER, stratum));
873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
889f6c5813SMatthew G. Knepley }
899f6c5813SMatthew G. Knepley
DMLabelSetUp_Ephemeral(DMLabel label)909f6c5813SMatthew G. Knepley static PetscErrorCode DMLabelSetUp_Ephemeral(DMLabel label)
919f6c5813SMatthew G. Knepley {
929f6c5813SMatthew G. Knepley DMLabel olabel;
939f6c5813SMatthew G. Knepley IS valueIS;
949f6c5813SMatthew G. Knepley const PetscInt *values;
959f6c5813SMatthew G. Knepley PetscInt defValue, Nv;
969f6c5813SMatthew G. Knepley
979f6c5813SMatthew G. Knepley PetscFunctionBegin;
989f6c5813SMatthew G. Knepley PetscCall(DMLabelEphemeralGetLabel(label, &olabel));
999f6c5813SMatthew G. Knepley PetscCall(DMLabelGetDefaultValue(olabel, &defValue));
1009f6c5813SMatthew G. Knepley PetscCall(DMLabelSetDefaultValue(label, defValue));
1019f6c5813SMatthew G. Knepley PetscCall(DMLabelGetNumValues(olabel, &Nv));
1029f6c5813SMatthew G. Knepley PetscCall(DMLabelGetValueIS(olabel, &valueIS));
1039f6c5813SMatthew G. Knepley PetscCall(ISGetIndices(valueIS, &values));
1049f6c5813SMatthew G. Knepley PetscCall(DMLabelAddStrataIS(label, valueIS));
1059f6c5813SMatthew G. Knepley for (PetscInt v = 0; v < Nv; ++v) PetscCall(DMLabelEphemeralComputeStratumSize_Private(label, values[v]));
1069f6c5813SMatthew G. Knepley PetscCall(ISRestoreIndices(valueIS, &values));
1079f6c5813SMatthew G. Knepley PetscCall(ISDestroy(&valueIS));
1089f6c5813SMatthew G. Knepley label->readonly = PETSC_TRUE;
1093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1109f6c5813SMatthew G. Knepley }
1119f6c5813SMatthew G. Knepley
DMLabelView_Ephemeral_Ascii(DMLabel label,PetscViewer viewer)1129f6c5813SMatthew G. Knepley static PetscErrorCode DMLabelView_Ephemeral_Ascii(DMLabel label, PetscViewer viewer)
1139f6c5813SMatthew G. Knepley {
1149f6c5813SMatthew G. Knepley DMLabel olabel;
1159f6c5813SMatthew G. Knepley PetscMPIInt rank;
1169f6c5813SMatthew G. Knepley
1179f6c5813SMatthew G. Knepley PetscFunctionBegin;
1189f6c5813SMatthew G. Knepley PetscCall(DMLabelEphemeralGetLabel(label, &olabel));
1199f6c5813SMatthew G. Knepley PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
1209f6c5813SMatthew G. Knepley PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1219f6c5813SMatthew G. Knepley if (olabel) {
1229f6c5813SMatthew G. Knepley IS valueIS;
1239f6c5813SMatthew G. Knepley const PetscInt *values;
1249f6c5813SMatthew G. Knepley PetscInt Nv, v;
1259f6c5813SMatthew G. Knepley const char *name;
1269f6c5813SMatthew G. Knepley
1279f6c5813SMatthew G. Knepley PetscCall(PetscObjectGetName((PetscObject)label, &name));
1289f6c5813SMatthew G. Knepley PetscCall(PetscViewerASCIIPrintf(viewer, "Ephemeral Label '%s':\n", name));
1299f6c5813SMatthew G. Knepley PetscCall(DMLabelGetNumValues(olabel, &Nv));
1309f6c5813SMatthew G. Knepley PetscCall(DMLabelGetValueIS(olabel, &valueIS));
1319f6c5813SMatthew G. Knepley PetscCall(ISGetIndices(valueIS, &values));
1329f6c5813SMatthew G. Knepley for (v = 0; v < Nv; ++v) {
1339f6c5813SMatthew G. Knepley IS pointIS;
1349f6c5813SMatthew G. Knepley const PetscInt value = values[v];
1359f6c5813SMatthew G. Knepley const PetscInt *points;
1369f6c5813SMatthew G. Knepley PetscInt n, p;
1379f6c5813SMatthew G. Knepley
1389f6c5813SMatthew G. Knepley PetscCall(DMLabelGetStratumIS(olabel, value, &pointIS));
1399f6c5813SMatthew G. Knepley PetscCall(ISGetIndices(pointIS, &points));
1409f6c5813SMatthew G. Knepley PetscCall(ISGetLocalSize(pointIS, &n));
1419f6c5813SMatthew G. Knepley for (p = 0; p < n; ++p) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %" PetscInt_FMT " (%" PetscInt_FMT ")\n", rank, points[p], value));
1429f6c5813SMatthew G. Knepley PetscCall(ISRestoreIndices(pointIS, &points));
1439f6c5813SMatthew G. Knepley PetscCall(ISDestroy(&pointIS));
1449f6c5813SMatthew G. Knepley }
1459f6c5813SMatthew G. Knepley PetscCall(ISRestoreIndices(valueIS, &values));
1469f6c5813SMatthew G. Knepley PetscCall(ISDestroy(&valueIS));
1479f6c5813SMatthew G. Knepley }
1489f6c5813SMatthew G. Knepley PetscCall(PetscViewerFlush(viewer));
1499f6c5813SMatthew G. Knepley PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1519f6c5813SMatthew G. Knepley }
1529f6c5813SMatthew G. Knepley
DMLabelView_Ephemeral(DMLabel label,PetscViewer viewer)1539f6c5813SMatthew G. Knepley static PetscErrorCode DMLabelView_Ephemeral(DMLabel label, PetscViewer viewer)
1549f6c5813SMatthew G. Knepley {
155*9f196a02SMartin Diehl PetscBool isascii;
1569f6c5813SMatthew G. Knepley
1579f6c5813SMatthew G. Knepley PetscFunctionBegin;
158*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
159*9f196a02SMartin Diehl if (isascii) PetscCall(DMLabelView_Ephemeral_Ascii(label, viewer));
1603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1619f6c5813SMatthew G. Knepley }
1629f6c5813SMatthew G. Knepley
DMLabelDuplicate_Ephemeral(DMLabel label,DMLabel * labelnew)1639f6c5813SMatthew G. Knepley static PetscErrorCode DMLabelDuplicate_Ephemeral(DMLabel label, DMLabel *labelnew)
1649f6c5813SMatthew G. Knepley {
1658dcf10e8SMatthew G. Knepley PetscObject obj;
1669f6c5813SMatthew G. Knepley
1679f6c5813SMatthew G. Knepley PetscFunctionBegin;
1688dcf10e8SMatthew G. Knepley PetscCall(PetscObjectQuery((PetscObject)label, "__original_label__", &obj));
1698dcf10e8SMatthew G. Knepley PetscCall(PetscObjectCompose((PetscObject)*labelnew, "__original_label__", obj));
1708dcf10e8SMatthew G. Knepley PetscCall(PetscObjectQuery((PetscObject)label, "__transform__", &obj));
1718dcf10e8SMatthew G. Knepley PetscCall(PetscObjectCompose((PetscObject)*labelnew, "__transform__", obj));
1728dcf10e8SMatthew G. Knepley PetscCall(DMLabelInitialize_Ephemeral(*labelnew));
1733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1749f6c5813SMatthew G. Knepley }
1759f6c5813SMatthew G. Knepley
DMLabelInitialize_Ephemeral(DMLabel label)17666976f2fSJacob Faibussowitsch static PetscErrorCode DMLabelInitialize_Ephemeral(DMLabel label)
1779f6c5813SMatthew G. Knepley {
1789f6c5813SMatthew G. Knepley PetscFunctionBegin;
1799f6c5813SMatthew G. Knepley label->ops->view = DMLabelView_Ephemeral;
1809f6c5813SMatthew G. Knepley label->ops->setup = DMLabelSetUp_Ephemeral;
1819f6c5813SMatthew G. Knepley label->ops->duplicate = DMLabelDuplicate_Ephemeral;
1829f6c5813SMatthew G. Knepley label->ops->getstratumis = DMLabelGetStratumIS_Ephemeral;
1833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1849f6c5813SMatthew G. Knepley }
1859f6c5813SMatthew G. Knepley
1869f6c5813SMatthew G. Knepley /*MC
1879f6c5813SMatthew G. Knepley DMLABELEPHEMERAL = "ephemeral" - This type of `DMLabel` is used to label ephemeral meshes.
1889f6c5813SMatthew G. Knepley
1899f6c5813SMatthew G. Knepley Ephemeral meshes are never concretely instantiated, but rather the answers to queries are created on the fly from a base mesh and a `DMPlexTransform` object.
1909f6c5813SMatthew G. Knepley For example, we could integrate over a refined mesh without ever storing the entire thing, just producing each cell closure one at a time. The label consists
1919f6c5813SMatthew G. Knepley of a base label and the same transform. Stratum are produced on demand, so that the time is in $O(max(M, N))$ where $M$ is the size of the original stratum,
1929f6c5813SMatthew G. Knepley and $N$ is the size of the output stratum. Queries take the same time, since we cannot sort the output.
1939f6c5813SMatthew G. Knepley
1949f6c5813SMatthew G. Knepley Level: intermediate
1959f6c5813SMatthew G. Knepley
1969f6c5813SMatthew G. Knepley .seealso: `DMLabel`, `DM`, `DMLabelType`, `DMLabelCreate()`, `DMLabelSetType()`
1979f6c5813SMatthew G. Knepley M*/
1989f6c5813SMatthew G. Knepley
DMLabelCreate_Ephemeral(DMLabel label)1999f6c5813SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLabelCreate_Ephemeral(DMLabel label)
2009f6c5813SMatthew G. Knepley {
2019f6c5813SMatthew G. Knepley PetscFunctionBegin;
2029f6c5813SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
2039f6c5813SMatthew G. Knepley PetscCall(DMLabelInitialize_Ephemeral(label));
2043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2059f6c5813SMatthew G. Knepley }
2069f6c5813SMatthew G. Knepley
2079f6c5813SMatthew G. Knepley /*@
2089f6c5813SMatthew G. Knepley DMLabelEphemeralGetLabel - Get the base label for this ephemeral label
2099f6c5813SMatthew G. Knepley
21020f4b53cSBarry Smith Not Collective
2119f6c5813SMatthew G. Knepley
2129f6c5813SMatthew G. Knepley Input Parameter:
21320f4b53cSBarry Smith . label - the `DMLabel`
2149f6c5813SMatthew G. Knepley
215aaa8cc7dSPierre Jolivet Output Parameter:
2169f6c5813SMatthew G. Knepley . olabel - the base label for this ephemeral label
2179f6c5813SMatthew G. Knepley
21820f4b53cSBarry Smith Level: intermediate
21920f4b53cSBarry Smith
2209f6c5813SMatthew G. Knepley Note:
2219f6c5813SMatthew G. Knepley Ephemeral labels are produced automatically from a base label and a `DMPlexTransform`.
2229f6c5813SMatthew G. Knepley
2239f6c5813SMatthew G. Knepley .seealso: `DMLabelEphemeralSetLabel()`, `DMLabelEphemeralGetTransform()`, `DMLabelSetType()`
2249f6c5813SMatthew G. Knepley @*/
DMLabelEphemeralGetLabel(DMLabel label,DMLabel * olabel)2259f6c5813SMatthew G. Knepley PetscErrorCode DMLabelEphemeralGetLabel(DMLabel label, DMLabel *olabel)
2269f6c5813SMatthew G. Knepley {
2279f6c5813SMatthew G. Knepley PetscFunctionBegin;
2289f6c5813SMatthew G. Knepley PetscCall(PetscObjectQuery((PetscObject)label, "__original_label__", (PetscObject *)olabel));
2293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2309f6c5813SMatthew G. Knepley }
2319f6c5813SMatthew G. Knepley
2329f6c5813SMatthew G. Knepley /*@
2339f6c5813SMatthew G. Knepley DMLabelEphemeralSetLabel - Set the base label for this ephemeral label
2349f6c5813SMatthew G. Knepley
23520f4b53cSBarry Smith Not Collective
2369f6c5813SMatthew G. Knepley
2379f6c5813SMatthew G. Knepley Input Parameters:
23820f4b53cSBarry Smith + label - the `DMLabel`
2399f6c5813SMatthew G. Knepley - olabel - the base label for this ephemeral label
2409f6c5813SMatthew G. Knepley
24120f4b53cSBarry Smith Level: intermediate
24220f4b53cSBarry Smith
2439f6c5813SMatthew G. Knepley Note:
2449f6c5813SMatthew G. Knepley Ephemeral labels are produced automatically from a base label and a `DMPlexTransform`.
2459f6c5813SMatthew G. Knepley
2469f6c5813SMatthew G. Knepley .seealso: `DMLabelEphemeralGetLabel()`, `DMLabelEphemeralSetTransform()`, `DMLabelSetType()`
2479f6c5813SMatthew G. Knepley @*/
DMLabelEphemeralSetLabel(DMLabel label,DMLabel olabel)2489f6c5813SMatthew G. Knepley PetscErrorCode DMLabelEphemeralSetLabel(DMLabel label, DMLabel olabel)
2499f6c5813SMatthew G. Knepley {
2509f6c5813SMatthew G. Knepley PetscFunctionBegin;
2519f6c5813SMatthew G. Knepley PetscCall(PetscObjectCompose((PetscObject)label, "__original_label__", (PetscObject)olabel));
2523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2539f6c5813SMatthew G. Knepley }
254