1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2*e8f14785SLisandro Dalcin #include <petsc/private/hashmapi.h> 3*e8f14785SLisandro Dalcin #include <petsc/private/hashmapij.h> 4*e8f14785SLisandro Dalcin 5*e8f14785SLisandro Dalcin /* HashIJKL */ 6*e8f14785SLisandro Dalcin 7*e8f14785SLisandro Dalcin #include <petsc/private/hashmap.h> 8*e8f14785SLisandro Dalcin 9*e8f14785SLisandro Dalcin typedef struct _PetscHashIJKLKey { PetscInt i, j, k, l; } PetscHashIJKLKey; 10*e8f14785SLisandro Dalcin 11*e8f14785SLisandro Dalcin #define PetscHashIJKLKeyHash(key) \ 12*e8f14785SLisandro Dalcin PetscHashCombine(PetscHashCombine(PetscHashInt((key).i),PetscHashInt((key).j)), \ 13*e8f14785SLisandro Dalcin PetscHashCombine(PetscHashInt((key).k),PetscHashInt((key).l))) 14*e8f14785SLisandro Dalcin 15*e8f14785SLisandro Dalcin #define PetscHashIJKLKeyEqual(k1,k2) \ 16*e8f14785SLisandro Dalcin (((k1).i==(k2).i) ? ((k1).j==(k2).j) ? ((k1).k==(k2).k) ? ((k1).l==(k2).l) : 0 : 0 : 0) 17*e8f14785SLisandro Dalcin 18*e8f14785SLisandro Dalcin PETSC_HASH_MAP(HashIJKL, PetscHashIJKLKey, PetscInt, PetscHashIJKLKeyHash, PetscHashIJKLKeyEqual, -1) 19*e8f14785SLisandro Dalcin 209f074e33SMatthew G Knepley 219f074e33SMatthew G Knepley /* 229a5b13a2SMatthew G Knepley DMPlexGetFaces_Internal - Gets groups of vertices that correspond to faces for the given cell 239f074e33SMatthew G Knepley */ 24439ece16SMatthew G. Knepley PetscErrorCode DMPlexGetFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 259f074e33SMatthew G Knepley { 269f074e33SMatthew G Knepley const PetscInt *cone = NULL; 279a5b13a2SMatthew G Knepley PetscInt maxConeSize, maxSupportSize, coneSize; 289f074e33SMatthew G Knepley PetscErrorCode ierr; 299f074e33SMatthew G Knepley 309f074e33SMatthew G Knepley PetscFunctionBegin; 319f074e33SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 329a5b13a2SMatthew G Knepley ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 339f074e33SMatthew G Knepley ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 349f074e33SMatthew G Knepley ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 35439ece16SMatthew G. Knepley ierr = DMPlexGetRawFaces_Internal(dm, dim, coneSize, cone, numFaces, faceSize, faces);CHKERRQ(ierr); 36439ece16SMatthew G. Knepley PetscFunctionReturn(0); 37439ece16SMatthew G. Knepley } 38439ece16SMatthew G. Knepley 39439ece16SMatthew G. Knepley /* 40439ece16SMatthew G. Knepley DMPlexRestoreFaces_Internal - Restores the array 41439ece16SMatthew G. Knepley */ 42439ece16SMatthew G. Knepley PetscErrorCode DMPlexRestoreFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 43439ece16SMatthew G. Knepley { 44439ece16SMatthew G. Knepley PetscErrorCode ierr; 45439ece16SMatthew G. Knepley 46439ece16SMatthew G. Knepley PetscFunctionBegin; 4769291d52SBarry Smith ierr = DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faces);CHKERRQ(ierr); 48439ece16SMatthew G. Knepley PetscFunctionReturn(0); 49439ece16SMatthew G. Knepley } 50439ece16SMatthew G. Knepley 51439ece16SMatthew G. Knepley /* 52439ece16SMatthew G. Knepley DMPlexGetRawFaces_Internal - Gets groups of vertices that correspond to faces for the given cone 53439ece16SMatthew G. Knepley */ 54439ece16SMatthew G. Knepley PetscErrorCode DMPlexGetRawFaces_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 55439ece16SMatthew G. Knepley { 56439ece16SMatthew G. Knepley PetscInt *facesTmp; 57439ece16SMatthew G. Knepley PetscInt maxConeSize, maxSupportSize; 58439ece16SMatthew G. Knepley PetscErrorCode ierr; 59439ece16SMatthew G. Knepley 60439ece16SMatthew G. Knepley PetscFunctionBegin; 61439ece16SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 62439ece16SMatthew G. Knepley ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 6369291d52SBarry Smith if (faces) {ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), MPIU_INT, &facesTmp);CHKERRQ(ierr);} 649f074e33SMatthew G Knepley switch (dim) { 65c49d9212SMatthew G. Knepley case 1: 66c49d9212SMatthew G. Knepley switch (coneSize) { 67c49d9212SMatthew G. Knepley case 2: 68c49d9212SMatthew G. Knepley if (faces) { 69c49d9212SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 70c49d9212SMatthew G. Knepley *faces = facesTmp; 71c49d9212SMatthew G. Knepley } 72c49d9212SMatthew G. Knepley if (numFaces) *numFaces = 2; 73c49d9212SMatthew G. Knepley if (faceSize) *faceSize = 1; 74c49d9212SMatthew G. Knepley break; 75c49d9212SMatthew G. Knepley default: 76c49d9212SMatthew G. Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 77c49d9212SMatthew G. Knepley } 78c49d9212SMatthew G. Knepley break; 799f074e33SMatthew G Knepley case 2: 809f074e33SMatthew G Knepley switch (coneSize) { 819f074e33SMatthew G Knepley case 3: 829a5b13a2SMatthew G Knepley if (faces) { 839a5b13a2SMatthew G Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 849a5b13a2SMatthew G Knepley facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 859a5b13a2SMatthew G Knepley facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 869a5b13a2SMatthew G Knepley *faces = facesTmp; 879a5b13a2SMatthew G Knepley } 889a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 3; 899a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 2; 909f074e33SMatthew G Knepley break; 919f074e33SMatthew G Knepley case 4: 929a5b13a2SMatthew G Knepley /* Vertices follow right hand rule */ 939a5b13a2SMatthew G Knepley if (faces) { 949a5b13a2SMatthew G Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 959a5b13a2SMatthew G Knepley facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 969a5b13a2SMatthew G Knepley facesTmp[4] = cone[2]; facesTmp[5] = cone[3]; 979a5b13a2SMatthew G Knepley facesTmp[6] = cone[3]; facesTmp[7] = cone[0]; 989a5b13a2SMatthew G Knepley *faces = facesTmp; 999a5b13a2SMatthew G Knepley } 1009a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 4; 1019a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 2; 1029f074e33SMatthew G Knepley break; 1039f074e33SMatthew G Knepley default: 1049f074e33SMatthew G Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 1059f074e33SMatthew G Knepley } 1069f074e33SMatthew G Knepley break; 1079f074e33SMatthew G Knepley case 3: 1089f074e33SMatthew G Knepley switch (coneSize) { 1099f074e33SMatthew G Knepley case 3: 1109a5b13a2SMatthew G Knepley if (faces) { 1119a5b13a2SMatthew G Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 1129a5b13a2SMatthew G Knepley facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 1139a5b13a2SMatthew G Knepley facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 1149a5b13a2SMatthew G Knepley *faces = facesTmp; 1159a5b13a2SMatthew G Knepley } 1169a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 3; 1179a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 2; 1189f074e33SMatthew G Knepley break; 1199f074e33SMatthew G Knepley case 4: 1202e1b13c2SMatthew G. Knepley /* Vertices of first face follow right hand rule and normal points away from last vertex */ 1219a5b13a2SMatthew G Knepley if (faces) { 1222e1b13c2SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; 1232e1b13c2SMatthew G. Knepley facesTmp[3] = cone[0]; facesTmp[4] = cone[3]; facesTmp[5] = cone[1]; 1242e1b13c2SMatthew G. Knepley facesTmp[6] = cone[0]; facesTmp[7] = cone[2]; facesTmp[8] = cone[3]; 1252e1b13c2SMatthew G. Knepley facesTmp[9] = cone[2]; facesTmp[10] = cone[1]; facesTmp[11] = cone[3]; 1269a5b13a2SMatthew G Knepley *faces = facesTmp; 1279a5b13a2SMatthew G Knepley } 1289a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 4; 1299a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 3; 1309a5b13a2SMatthew G Knepley break; 1319a5b13a2SMatthew G Knepley case 8: 132e368ef20SMatthew G. Knepley /* 7--------6 133e368ef20SMatthew G. Knepley /| /| 134e368ef20SMatthew G. Knepley / | / | 135e368ef20SMatthew G. Knepley 4--------5 | 136e368ef20SMatthew G. Knepley | | | | 137e368ef20SMatthew G. Knepley | | | | 138e368ef20SMatthew G. Knepley | 1--------2 139e368ef20SMatthew G. Knepley | / | / 140e368ef20SMatthew G. Knepley |/ |/ 141e368ef20SMatthew G. Knepley 0--------3 142e368ef20SMatthew G. Knepley */ 1439a5b13a2SMatthew G Knepley if (faces) { 14451cfd6a4SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; facesTmp[3] = cone[3]; /* Bottom */ 14551cfd6a4SMatthew G. Knepley facesTmp[4] = cone[4]; facesTmp[5] = cone[5]; facesTmp[6] = cone[6]; facesTmp[7] = cone[7]; /* Top */ 14651cfd6a4SMatthew G. Knepley facesTmp[8] = cone[0]; facesTmp[9] = cone[3]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4]; /* Front */ 14751cfd6a4SMatthew G. Knepley facesTmp[12] = cone[2]; facesTmp[13] = cone[1]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6]; /* Back */ 14851cfd6a4SMatthew G. Knepley facesTmp[16] = cone[3]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5]; /* Right */ 14951cfd6a4SMatthew G. Knepley facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[1]; /* Left */ 1509a5b13a2SMatthew G Knepley *faces = facesTmp; 1519a5b13a2SMatthew G Knepley } 1529a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 6; 1539a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 4; 1549f074e33SMatthew G Knepley break; 1559f074e33SMatthew G Knepley default: 1569f074e33SMatthew G Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 1579f074e33SMatthew G Knepley } 1589f074e33SMatthew G Knepley break; 1599f074e33SMatthew G Knepley default: 1609f074e33SMatthew G Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim); 1619f074e33SMatthew G Knepley } 1629f074e33SMatthew G Knepley PetscFunctionReturn(0); 1639f074e33SMatthew G Knepley } 1649f074e33SMatthew G Knepley 1659a5b13a2SMatthew G Knepley /* This interpolates faces for cells at some stratum */ 1669a5b13a2SMatthew G Knepley static PetscErrorCode DMPlexInterpolateFaces_Internal(DM dm, PetscInt cellDepth, DM idm) 1679f074e33SMatthew G Knepley { 168d4efc6f1SMatthew G. Knepley DMLabel subpointMap; 1699a5b13a2SMatthew G Knepley PetscHashIJKL faceTable; 1709a5b13a2SMatthew G Knepley PetscInt *pStart, *pEnd; 1719a5b13a2SMatthew G Knepley PetscInt cellDim, depth, faceDepth = cellDepth, numPoints = 0, faceSizeAll = 0, face, c, d; 1729f074e33SMatthew G Knepley PetscErrorCode ierr; 1739f074e33SMatthew G Knepley 1749f074e33SMatthew G Knepley PetscFunctionBegin; 175c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &cellDim);CHKERRQ(ierr); 176d4efc6f1SMatthew G. Knepley /* HACK: I need a better way to determine face dimension, or an alternative to GetFaces() */ 177d4efc6f1SMatthew G. Knepley ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr); 178d4efc6f1SMatthew G. Knepley if (subpointMap) ++cellDim; 1799a5b13a2SMatthew G Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1809a5b13a2SMatthew G Knepley ++depth; 1819a5b13a2SMatthew G Knepley ++cellDepth; 1829a5b13a2SMatthew G Knepley cellDim -= depth - cellDepth; 183dcca6d9dSJed Brown ierr = PetscMalloc2(depth+1,&pStart,depth+1,&pEnd);CHKERRQ(ierr); 1849a5b13a2SMatthew G Knepley for (d = depth-1; d >= faceDepth; --d) { 1859a5b13a2SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, d, &pStart[d+1], &pEnd[d+1]);CHKERRQ(ierr); 1869f074e33SMatthew G Knepley } 1879a5b13a2SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, -1, NULL, &pStart[faceDepth]);CHKERRQ(ierr); 1889a5b13a2SMatthew G Knepley pEnd[faceDepth] = pStart[faceDepth]; 1899a5b13a2SMatthew G Knepley for (d = faceDepth-1; d >= 0; --d) { 1909a5b13a2SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 1919f074e33SMatthew G Knepley } 1929a5b13a2SMatthew G Knepley if (pEnd[cellDepth] > pStart[cellDepth]) {ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr);} 1939a5b13a2SMatthew G Knepley if (faceSizeAll > 4) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll); 1949a5b13a2SMatthew G Knepley ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 1959a5b13a2SMatthew G Knepley for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) { 1969f074e33SMatthew G Knepley const PetscInt *cellFaces; 197735a0255SMatthew G. Knepley PetscInt numCellFaces, faceSize, cf; 1989f074e33SMatthew G Knepley 1999a5b13a2SMatthew G Knepley ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 2009a5b13a2SMatthew G Knepley if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 2019f074e33SMatthew G Knepley for (cf = 0; cf < numCellFaces; ++cf) { 2029a5b13a2SMatthew G Knepley const PetscInt *cellFace = &cellFaces[cf*faceSize]; 2039a5b13a2SMatthew G Knepley PetscHashIJKLKey key; 204*e8f14785SLisandro Dalcin PetscHashIter iter; 205*e8f14785SLisandro Dalcin PetscBool missing; 2069f074e33SMatthew G Knepley 2079a5b13a2SMatthew G Knepley if (faceSize == 2) { 2089a5b13a2SMatthew G Knepley key.i = PetscMin(cellFace[0], cellFace[1]); 2099a5b13a2SMatthew G Knepley key.j = PetscMax(cellFace[0], cellFace[1]); 2107fcd4ef0SJed Brown key.k = 0; 2117fcd4ef0SJed Brown key.l = 0; 2129a5b13a2SMatthew G Knepley } else { 2139a5b13a2SMatthew G Knepley key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0; 214302440fdSBarry Smith ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr); 2159f074e33SMatthew G Knepley } 216*e8f14785SLisandro Dalcin ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr); 217*e8f14785SLisandro Dalcin if (missing) {ierr = PetscHashIJKLIterSet(faceTable, iter, face++);CHKERRQ(ierr);} 2189a5b13a2SMatthew G Knepley } 219439ece16SMatthew G. Knepley ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 2209a5b13a2SMatthew G Knepley } 2219a5b13a2SMatthew G Knepley pEnd[faceDepth] = face; 2229a5b13a2SMatthew G Knepley ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 2239a5b13a2SMatthew G Knepley /* Count new points */ 2249a5b13a2SMatthew G Knepley for (d = 0; d <= depth; ++d) { 2259a5b13a2SMatthew G Knepley numPoints += pEnd[d]-pStart[d]; 2269a5b13a2SMatthew G Knepley } 2279a5b13a2SMatthew G Knepley ierr = DMPlexSetChart(idm, 0, numPoints);CHKERRQ(ierr); 2289a5b13a2SMatthew G Knepley /* Set cone sizes */ 2299a5b13a2SMatthew G Knepley for (d = 0; d <= depth; ++d) { 2309a5b13a2SMatthew G Knepley PetscInt coneSize, p; 2319f074e33SMatthew G Knepley 2329a5b13a2SMatthew G Knepley if (d == faceDepth) { 2339a5b13a2SMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 2349a5b13a2SMatthew G Knepley /* I see no way to do this if we admit faces of different shapes */ 2359a5b13a2SMatthew G Knepley ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr); 2369a5b13a2SMatthew G Knepley } 237a014044eSMatthew G Knepley } else if (d == cellDepth) { 238a014044eSMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 239a014044eSMatthew G Knepley /* Number of cell faces may be different from number of cell vertices*/ 240a014044eSMatthew G Knepley ierr = DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);CHKERRQ(ierr); 241a014044eSMatthew G Knepley ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 242a014044eSMatthew G Knepley } 2439a5b13a2SMatthew G Knepley } else { 2449a5b13a2SMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 2459a5b13a2SMatthew G Knepley ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 2469a5b13a2SMatthew G Knepley ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 2479f074e33SMatthew G Knepley } 2489f074e33SMatthew G Knepley } 2499f074e33SMatthew G Knepley } 2509f074e33SMatthew G Knepley ierr = DMSetUp(idm);CHKERRQ(ierr); 2519f074e33SMatthew G Knepley /* Get face cones from subsets of cell vertices */ 2529a5b13a2SMatthew G Knepley if (faceSizeAll > 4) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll); 2539a5b13a2SMatthew G Knepley ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 2549a5b13a2SMatthew G Knepley for (d = depth; d > cellDepth; --d) { 2559a5b13a2SMatthew G Knepley const PetscInt *cone; 2569a5b13a2SMatthew G Knepley PetscInt p; 2579a5b13a2SMatthew G Knepley 2589a5b13a2SMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 2599a5b13a2SMatthew G Knepley ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 2609a5b13a2SMatthew G Knepley ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr); 2619a5b13a2SMatthew G Knepley ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr); 2629a5b13a2SMatthew G Knepley ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr); 2639f074e33SMatthew G Knepley } 2649a5b13a2SMatthew G Knepley } 2659a5b13a2SMatthew G Knepley for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) { 2669f074e33SMatthew G Knepley const PetscInt *cellFaces; 267735a0255SMatthew G. Knepley PetscInt numCellFaces, faceSize, cf; 2689f074e33SMatthew G Knepley 2699a5b13a2SMatthew G Knepley ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 2709a5b13a2SMatthew G Knepley if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 2719f074e33SMatthew G Knepley for (cf = 0; cf < numCellFaces; ++cf) { 2729a5b13a2SMatthew G Knepley const PetscInt *cellFace = &cellFaces[cf*faceSize]; 2739a5b13a2SMatthew G Knepley PetscHashIJKLKey key; 274*e8f14785SLisandro Dalcin PetscHashIter iter; 275*e8f14785SLisandro Dalcin PetscBool missing; 2769f074e33SMatthew G Knepley 2779a5b13a2SMatthew G Knepley if (faceSize == 2) { 2789a5b13a2SMatthew G Knepley key.i = PetscMin(cellFace[0], cellFace[1]); 2799a5b13a2SMatthew G Knepley key.j = PetscMax(cellFace[0], cellFace[1]); 2807fcd4ef0SJed Brown key.k = 0; 2817fcd4ef0SJed Brown key.l = 0; 2829a5b13a2SMatthew G Knepley } else { 283*e8f14785SLisandro Dalcin key.i = cellFace[0]; 284*e8f14785SLisandro Dalcin key.j = cellFace[1]; 285*e8f14785SLisandro Dalcin key.k = cellFace[2]; 286*e8f14785SLisandro Dalcin key.l = faceSize > 3 ? cellFace[3] : 0; 287302440fdSBarry Smith ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr); 2889f074e33SMatthew G Knepley } 289*e8f14785SLisandro Dalcin ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr); 290735a0255SMatthew G. Knepley if (missing) { 2919a5b13a2SMatthew G Knepley ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr); 292*e8f14785SLisandro Dalcin ierr = PetscHashIJKLIterSet(faceTable, iter, face);CHKERRQ(ierr); 293735a0255SMatthew G. Knepley ierr = DMPlexInsertCone(idm, c, cf, face++);CHKERRQ(ierr); 2949a5b13a2SMatthew G Knepley } else { 2959a5b13a2SMatthew G Knepley const PetscInt *cone; 296735a0255SMatthew G. Knepley PetscInt coneSize, ornt, i, j, f; 2979f074e33SMatthew G Knepley 298*e8f14785SLisandro Dalcin ierr = PetscHashIJKLIterGet(faceTable, iter, &f);CHKERRQ(ierr); 2999a5b13a2SMatthew G Knepley ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr); 3002e1b13c2SMatthew G. Knepley /* Orient face: Do not allow reverse orientation at the first vertex */ 3019f074e33SMatthew G Knepley ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr); 3029f074e33SMatthew G Knepley ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr); 3039a5b13a2SMatthew G Knepley if (coneSize != faceSize) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of face vertices %D for face %D should be %D", coneSize, f, faceSize); 3049a5b13a2SMatthew G Knepley /* - First find the initial vertex */ 3059a5b13a2SMatthew G Knepley for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break; 3069a5b13a2SMatthew G Knepley /* - Try forward comparison */ 3079a5b13a2SMatthew G Knepley for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break; 3089a5b13a2SMatthew G Knepley if (j == faceSize) { 3099a5b13a2SMatthew G Knepley if ((faceSize == 2) && (i == 1)) ornt = -2; 3109a5b13a2SMatthew G Knepley else ornt = i; 3119a5b13a2SMatthew G Knepley } else { 3129a5b13a2SMatthew G Knepley /* - Try backward comparison */ 3139a5b13a2SMatthew G Knepley for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break; 3142e1b13c2SMatthew G. Knepley if (j == faceSize) { 3152e1b13c2SMatthew G. Knepley if (i == 0) ornt = -faceSize; 316dc1a705cSMatthew G. Knepley else ornt = -i; 3172e1b13c2SMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation"); 3189a5b13a2SMatthew G Knepley } 3199a5b13a2SMatthew G Knepley ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr); 3209f074e33SMatthew G Knepley } 3219f074e33SMatthew G Knepley } 322439ece16SMatthew G. Knepley ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 3239f074e33SMatthew G Knepley } 3249a5b13a2SMatthew G Knepley if (face != pEnd[faceDepth]) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Invalid number of faces %D should be %D", face-pStart[faceDepth], pEnd[faceDepth]-pStart[faceDepth]); 325c907b753SJed Brown ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 3269a5b13a2SMatthew G Knepley ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 3276551a8c7SMatthew G. Knepley ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 3289a5b13a2SMatthew G Knepley ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr); 3299a5b13a2SMatthew G Knepley ierr = DMPlexStratify(idm);CHKERRQ(ierr); 3309f074e33SMatthew G Knepley PetscFunctionReturn(0); 3319f074e33SMatthew G Knepley } 3329f074e33SMatthew G Knepley 3337bffde78SMichael Lange /* This interpolates the PointSF in parallel following local interpolation */ 3347bffde78SMichael Lange static PetscErrorCode DMPlexInterpolatePointSF(DM dm, PetscSF pointSF, PetscInt depth) 3357bffde78SMichael Lange { 3369852e123SBarry Smith PetscMPIInt size, rank; 3377bffde78SMichael Lange PetscInt p, c, d, dof, offset; 338ffd6914dSSatish Balay PetscInt numLeaves, numRoots, candidatesSize, candidatesRemoteSize; 3397bffde78SMichael Lange const PetscInt *localPoints; 3407bffde78SMichael Lange const PetscSFNode *remotePoints; 3417bffde78SMichael Lange PetscSFNode *candidates, *candidatesRemote, *claims; 3427bffde78SMichael Lange PetscSection candidateSection, candidateSectionRemote, claimSection; 343*e8f14785SLisandro Dalcin PetscHMapI leafhash; 344*e8f14785SLisandro Dalcin PetscHMapIJ roothash; 3457bffde78SMichael Lange PetscHashIJKey key; 3467bffde78SMichael Lange PetscErrorCode ierr; 3477bffde78SMichael Lange 3487bffde78SMichael Lange PetscFunctionBegin; 3499852e123SBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr); 3507bffde78SMichael Lange ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 3517bffde78SMichael Lange ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 3529852e123SBarry Smith if (size < 2 || numRoots < 0) PetscFunctionReturn(0); 35325afeb17SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr); 3547bffde78SMichael Lange /* Build hashes of points in the SF for efficient lookup */ 355*e8f14785SLisandro Dalcin ierr = PetscHMapICreate(&leafhash);CHKERRQ(ierr); 356*e8f14785SLisandro Dalcin ierr = PetscHMapIJCreate(&roothash);CHKERRQ(ierr); 3577bffde78SMichael Lange for (p = 0; p < numLeaves; ++p) { 358*e8f14785SLisandro Dalcin ierr = PetscHMapISet(leafhash, localPoints[p], p);CHKERRQ(ierr); 359*e8f14785SLisandro Dalcin key.i = remotePoints[p].index; 360*e8f14785SLisandro Dalcin key.j = remotePoints[p].rank; 361*e8f14785SLisandro Dalcin ierr = PetscHMapIJSet(roothash, key, p);CHKERRQ(ierr); 3627bffde78SMichael Lange } 3637bffde78SMichael Lange /* Build a section / SFNode array of candidate points in the single-level adjacency of leaves, 3647bffde78SMichael Lange where each candidate is defined by the root entry for the other vertex that defines the edge. */ 3657bffde78SMichael Lange ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSection);CHKERRQ(ierr); 3667bffde78SMichael Lange ierr = PetscSectionSetChart(candidateSection, 0, numRoots);CHKERRQ(ierr); 3677bffde78SMichael Lange { 368ffd6914dSSatish Balay PetscInt leaf, root, idx, a, *adj = NULL; 3697bffde78SMichael Lange for (p = 0; p < numLeaves; ++p) { 3707bffde78SMichael Lange PetscInt adjSize = PETSC_DETERMINE; 3717bffde78SMichael Lange ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr); 3727bffde78SMichael Lange for (a = 0; a < adjSize; ++a) { 373*e8f14785SLisandro Dalcin ierr = PetscHMapIGet(leafhash, adj[a], &leaf);CHKERRQ(ierr); 3747bffde78SMichael Lange if (leaf >= 0) {ierr = PetscSectionAddDof(candidateSection, localPoints[p], 1);CHKERRQ(ierr);} 3757bffde78SMichael Lange } 3767bffde78SMichael Lange } 3777bffde78SMichael Lange ierr = PetscSectionSetUp(candidateSection);CHKERRQ(ierr); 3787bffde78SMichael Lange ierr = PetscSectionGetStorageSize(candidateSection, &candidatesSize);CHKERRQ(ierr); 3797bffde78SMichael Lange ierr = PetscMalloc1(candidatesSize, &candidates);CHKERRQ(ierr); 3807bffde78SMichael Lange for (p = 0; p < numLeaves; ++p) { 3817bffde78SMichael Lange PetscInt adjSize = PETSC_DETERMINE; 3827bffde78SMichael Lange ierr = PetscSectionGetOffset(candidateSection, localPoints[p], &offset);CHKERRQ(ierr); 3837bffde78SMichael Lange ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr); 3847bffde78SMichael Lange for (idx = 0, a = 0; a < adjSize; ++a) { 385*e8f14785SLisandro Dalcin ierr = PetscHMapIGet(leafhash, adj[a], &root);CHKERRQ(ierr); 3867bffde78SMichael Lange if (root >= 0) candidates[offset+idx++] = remotePoints[root]; 3877bffde78SMichael Lange } 3887bffde78SMichael Lange } 3897bffde78SMichael Lange ierr = PetscFree(adj);CHKERRQ(ierr); 3907bffde78SMichael Lange } 3917bffde78SMichael Lange /* Gather candidate section / array pair into the root partition via inverse(multi(pointSF)). */ 3927bffde78SMichael Lange { 3937bffde78SMichael Lange PetscSF sfMulti, sfInverse, sfCandidates; 3947bffde78SMichael Lange PetscInt *remoteOffsets; 3957bffde78SMichael Lange ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr); 3967bffde78SMichael Lange ierr = PetscSFCreateInverseSF(sfMulti, &sfInverse);CHKERRQ(ierr); 3977bffde78SMichael Lange ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSectionRemote);CHKERRQ(ierr); 3987bffde78SMichael Lange ierr = PetscSFDistributeSection(sfInverse, candidateSection, &remoteOffsets, candidateSectionRemote);CHKERRQ(ierr); 3997bffde78SMichael Lange ierr = PetscSFCreateSectionSF(sfInverse, candidateSection, remoteOffsets, candidateSectionRemote, &sfCandidates);CHKERRQ(ierr); 4007bffde78SMichael Lange ierr = PetscSectionGetStorageSize(candidateSectionRemote, &candidatesRemoteSize);CHKERRQ(ierr); 4017bffde78SMichael Lange ierr = PetscMalloc1(candidatesRemoteSize, &candidatesRemote);CHKERRQ(ierr); 4027bffde78SMichael Lange ierr = PetscSFBcastBegin(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr); 4037bffde78SMichael Lange ierr = PetscSFBcastEnd(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr); 4047bffde78SMichael Lange ierr = PetscSFDestroy(&sfInverse);CHKERRQ(ierr); 4057bffde78SMichael Lange ierr = PetscSFDestroy(&sfCandidates);CHKERRQ(ierr); 4067bffde78SMichael Lange ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 4077bffde78SMichael Lange } 4087bffde78SMichael Lange /* Walk local roots and check for each remote candidate whether we know all required points, 4097bffde78SMichael Lange either from owning it or having a root entry in the point SF. If we do we place a claim 4107bffde78SMichael Lange by replacing the vertex number with our edge ID. */ 4117bffde78SMichael Lange { 4127bffde78SMichael Lange PetscInt idx, root, joinSize, vertices[2]; 4137bffde78SMichael Lange const PetscInt *rootdegree, *join = NULL; 4147bffde78SMichael Lange ierr = PetscSFComputeDegreeBegin(pointSF, &rootdegree);CHKERRQ(ierr); 4157bffde78SMichael Lange ierr = PetscSFComputeDegreeEnd(pointSF, &rootdegree);CHKERRQ(ierr); 4167bffde78SMichael Lange /* Loop remote edge connections and put in a claim if both vertices are known */ 4177bffde78SMichael Lange for (idx = 0, p = 0; p < numRoots; ++p) { 4187bffde78SMichael Lange for (d = 0; d < rootdegree[p]; ++d) { 4197bffde78SMichael Lange ierr = PetscSectionGetDof(candidateSectionRemote, idx, &dof);CHKERRQ(ierr); 4207bffde78SMichael Lange ierr = PetscSectionGetOffset(candidateSectionRemote, idx, &offset);CHKERRQ(ierr); 4217bffde78SMichael Lange for (c = 0; c < dof; ++c) { 4227bffde78SMichael Lange /* We own both vertices, so we claim the edge by replacing vertex with edge */ 4237bffde78SMichael Lange if (candidatesRemote[offset+c].rank == rank) { 4247bffde78SMichael Lange vertices[0] = p; vertices[1] = candidatesRemote[offset+c].index; 4257bffde78SMichael Lange ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 4267bffde78SMichael Lange if (joinSize == 1) candidatesRemote[offset+c].index = join[0]; 4277bffde78SMichael Lange ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 4287bffde78SMichael Lange continue; 4297bffde78SMichael Lange } 4307bffde78SMichael Lange /* If we own one vertex and share a root with the other, we claim it */ 431*e8f14785SLisandro Dalcin key.i = candidatesRemote[offset+c].index; 432*e8f14785SLisandro Dalcin key.j = candidatesRemote[offset+c].rank; 433*e8f14785SLisandro Dalcin ierr = PetscHMapIJGet(roothash, key, &root);CHKERRQ(ierr); 4347bffde78SMichael Lange if (root >= 0) { 4357bffde78SMichael Lange vertices[0] = p; vertices[1] = localPoints[root]; 4367bffde78SMichael Lange ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 4377bffde78SMichael Lange if (joinSize == 1) { 4387bffde78SMichael Lange candidatesRemote[offset+c].index = join[0]; 4397bffde78SMichael Lange candidatesRemote[offset+c].rank = rank; 4407bffde78SMichael Lange } 4417bffde78SMichael Lange ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 4427bffde78SMichael Lange } 4437bffde78SMichael Lange } 4447bffde78SMichael Lange idx++; 4457bffde78SMichael Lange } 4467bffde78SMichael Lange } 4477bffde78SMichael Lange } 4487bffde78SMichael Lange /* Push claims back to receiver via the MultiSF and derive new pointSF mapping on receiver */ 4497bffde78SMichael Lange { 4507bffde78SMichael Lange PetscSF sfMulti, sfClaims, sfPointNew; 451*e8f14785SLisandro Dalcin PetscHMapI claimshash; 4527bffde78SMichael Lange PetscInt size, pStart, pEnd, root, joinSize, numLocalNew; 4537bffde78SMichael Lange PetscInt *remoteOffsets, *localPointsNew, vertices[2]; 454ffd6914dSSatish Balay const PetscInt *join = NULL; 4557bffde78SMichael Lange PetscSFNode *remotePointsNew; 4567bffde78SMichael Lange ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr); 4577bffde78SMichael Lange ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &claimSection);CHKERRQ(ierr); 4587bffde78SMichael Lange ierr = PetscSFDistributeSection(sfMulti, candidateSectionRemote, &remoteOffsets, claimSection);CHKERRQ(ierr); 4597bffde78SMichael Lange ierr = PetscSFCreateSectionSF(sfMulti, candidateSectionRemote, remoteOffsets, claimSection, &sfClaims);CHKERRQ(ierr); 4607bffde78SMichael Lange ierr = PetscSectionGetStorageSize(claimSection, &size);CHKERRQ(ierr); 4617bffde78SMichael Lange ierr = PetscMalloc1(size, &claims);CHKERRQ(ierr); 4627bffde78SMichael Lange ierr = PetscSFBcastBegin(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr); 4637bffde78SMichael Lange ierr = PetscSFBcastEnd(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr); 4647bffde78SMichael Lange ierr = PetscSFDestroy(&sfClaims);CHKERRQ(ierr); 4657bffde78SMichael Lange ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 4667bffde78SMichael Lange /* Walk the original section of local supports and add an SF entry for each updated item */ 467*e8f14785SLisandro Dalcin ierr = PetscHMapICreate(&claimshash);CHKERRQ(ierr); 4687bffde78SMichael Lange for (p = 0; p < numRoots; ++p) { 4697bffde78SMichael Lange ierr = PetscSectionGetDof(candidateSection, p, &dof);CHKERRQ(ierr); 4707bffde78SMichael Lange ierr = PetscSectionGetOffset(candidateSection, p, &offset);CHKERRQ(ierr); 4717bffde78SMichael Lange for (d = 0; d < dof; ++d) { 4727bffde78SMichael Lange if (candidates[offset+d].index != claims[offset+d].index) { 473*e8f14785SLisandro Dalcin key.i = candidates[offset+d].index; 474*e8f14785SLisandro Dalcin key.j = candidates[offset+d].rank; 475*e8f14785SLisandro Dalcin ierr = PetscHMapIJGet(roothash, key, &root);CHKERRQ(ierr); 4767bffde78SMichael Lange if (root >= 0) { 4777bffde78SMichael Lange vertices[0] = p; vertices[1] = localPoints[root]; 4787bffde78SMichael Lange ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 479*e8f14785SLisandro Dalcin if (joinSize == 1) {ierr = PetscHMapISet(claimshash, join[0], offset+d);CHKERRQ(ierr);} 4807bffde78SMichael Lange ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 4817bffde78SMichael Lange } 4827bffde78SMichael Lange } 4837bffde78SMichael Lange } 4847bffde78SMichael Lange } 4857bffde78SMichael Lange /* Create new pointSF from hashed claims */ 486*e8f14785SLisandro Dalcin ierr = PetscHMapIGetSize(claimshash, &numLocalNew);CHKERRQ(ierr); 4877bffde78SMichael Lange ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 4887bffde78SMichael Lange ierr = PetscMalloc1(numLeaves + numLocalNew, &localPointsNew);CHKERRQ(ierr); 4897bffde78SMichael Lange ierr = PetscMalloc1(numLeaves + numLocalNew, &remotePointsNew);CHKERRQ(ierr); 4907bffde78SMichael Lange for (p = 0; p < numLeaves; ++p) { 4917bffde78SMichael Lange localPointsNew[p] = localPoints[p]; 4927bffde78SMichael Lange remotePointsNew[p].index = remotePoints[p].index; 4937bffde78SMichael Lange remotePointsNew[p].rank = remotePoints[p].rank; 4947bffde78SMichael Lange } 495f3190fc2SToby Isaac p = numLeaves; 496*e8f14785SLisandro Dalcin ierr = PetscHMapIGetKeys(claimshash, &p, localPointsNew);CHKERRQ(ierr); 497f3190fc2SToby Isaac ierr = PetscSortInt(numLocalNew, &localPointsNew[numLeaves]);CHKERRQ(ierr); 4987bffde78SMichael Lange for (p = numLeaves; p < numLeaves + numLocalNew; ++p) { 499*e8f14785SLisandro Dalcin ierr = PetscHMapIGet(claimshash, localPointsNew[p], &offset);CHKERRQ(ierr); 5007bffde78SMichael Lange remotePointsNew[p] = claims[offset]; 5017bffde78SMichael Lange } 5027bffde78SMichael Lange ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), &sfPointNew);CHKERRQ(ierr); 5037bffde78SMichael Lange ierr = PetscSFSetGraph(sfPointNew, pEnd-pStart, numLeaves+numLocalNew, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr); 5047bffde78SMichael Lange ierr = DMSetPointSF(dm, sfPointNew);CHKERRQ(ierr); 5057bffde78SMichael Lange ierr = PetscSFDestroy(&sfPointNew);CHKERRQ(ierr); 506*e8f14785SLisandro Dalcin ierr = PetscHMapIDestroy(&claimshash);CHKERRQ(ierr); 5077bffde78SMichael Lange } 508*e8f14785SLisandro Dalcin ierr = PetscHMapIDestroy(&leafhash);CHKERRQ(ierr); 509*e8f14785SLisandro Dalcin ierr = PetscHMapIJDestroy(&roothash);CHKERRQ(ierr); 5107bffde78SMichael Lange ierr = PetscSectionDestroy(&candidateSection);CHKERRQ(ierr); 5117bffde78SMichael Lange ierr = PetscSectionDestroy(&candidateSectionRemote);CHKERRQ(ierr); 5127bffde78SMichael Lange ierr = PetscSectionDestroy(&claimSection);CHKERRQ(ierr); 5137bffde78SMichael Lange ierr = PetscFree(candidates);CHKERRQ(ierr); 5147bffde78SMichael Lange ierr = PetscFree(candidatesRemote);CHKERRQ(ierr); 5157bffde78SMichael Lange ierr = PetscFree(claims);CHKERRQ(ierr); 51625afeb17SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr); 5177bffde78SMichael Lange PetscFunctionReturn(0); 5187bffde78SMichael Lange } 5197bffde78SMichael Lange 5200c795ddcSMatthew G. Knepley /*@C 52180330477SMatthew G. Knepley DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc. 52280330477SMatthew G. Knepley 52380330477SMatthew G. Knepley Collective on DM 52480330477SMatthew G. Knepley 525e56d480eSMatthew G. Knepley Input Parameters: 526e56d480eSMatthew G. Knepley + dm - The DMPlex object with only cells and vertices 52710a67516SMatthew G. Knepley - dmInt - The interpolated DM 52880330477SMatthew G. Knepley 52980330477SMatthew G. Knepley Output Parameter: 5304e3744c5SMatthew G. Knepley . dmInt - The complete DMPlex object 53180330477SMatthew G. Knepley 53280330477SMatthew G. Knepley Level: intermediate 53380330477SMatthew G. Knepley 53495452b02SPatrick Sanan Notes: 53595452b02SPatrick Sanan It does not copy over the coordinates. 53643eeeb2dSStefano Zampini 53780330477SMatthew G. Knepley .keywords: mesh 53843eeeb2dSStefano Zampini .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates() 53980330477SMatthew G. Knepley @*/ 5409f074e33SMatthew G Knepley PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt) 5419f074e33SMatthew G Knepley { 5429a5b13a2SMatthew G Knepley DM idm, odm = dm; 5437bffde78SMichael Lange PetscSF sfPoint; 5442e1b13c2SMatthew G. Knepley PetscInt depth, dim, d; 54510a67516SMatthew G. Knepley const char *name; 5469f074e33SMatthew G Knepley PetscErrorCode ierr; 5479f074e33SMatthew G Knepley 5489f074e33SMatthew G Knepley PetscFunctionBegin; 54910a67516SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 55010a67516SMatthew G. Knepley PetscValidPointer(dmInt, 2); 551a72f3261SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr); 5522e1b13c2SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 553c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 554b21b8912SMatthew G. Knepley if ((depth == dim) || (dim <= 1)) { 55576b791aaSMatthew G. Knepley ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 55676b791aaSMatthew G. Knepley idm = dm; 557b21b8912SMatthew G. Knepley } else { 5589a5b13a2SMatthew G Knepley for (d = 1; d < dim; ++d) { 5599a5b13a2SMatthew G Knepley /* Create interpolated mesh */ 56010a67516SMatthew G. Knepley ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr); 5619a5b13a2SMatthew G Knepley ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr); 562c73cfb54SMatthew G. Knepley ierr = DMSetDimension(idm, dim);CHKERRQ(ierr); 5637bffde78SMichael Lange if (depth > 0) { 5647bffde78SMichael Lange ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr); 5657bffde78SMichael Lange ierr = DMGetPointSF(odm, &sfPoint);CHKERRQ(ierr); 5667bffde78SMichael Lange ierr = DMPlexInterpolatePointSF(idm, sfPoint, depth);CHKERRQ(ierr); 5677bffde78SMichael Lange } 5689a5b13a2SMatthew G Knepley if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);} 5699a5b13a2SMatthew G Knepley odm = idm; 5709f074e33SMatthew G Knepley } 57110a67516SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) dm, &name);CHKERRQ(ierr); 57210a67516SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) idm, name);CHKERRQ(ierr); 57310a67516SMatthew G. Knepley ierr = DMPlexCopyCoordinates(dm, idm);CHKERRQ(ierr); 57410a67516SMatthew G. Knepley ierr = DMCopyLabels(dm, idm);CHKERRQ(ierr); 57584699f70SSatish Balay } 57643eeeb2dSStefano Zampini { 57743eeeb2dSStefano Zampini PetscBool isper; 57843eeeb2dSStefano Zampini const PetscReal *maxCell, *L; 57943eeeb2dSStefano Zampini const DMBoundaryType *bd; 58043eeeb2dSStefano Zampini 58143eeeb2dSStefano Zampini ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr); 58243eeeb2dSStefano Zampini ierr = DMSetPeriodicity(idm,isper,maxCell,L,bd);CHKERRQ(ierr); 58343eeeb2dSStefano Zampini } 5849a5b13a2SMatthew G Knepley *dmInt = idm; 585a72f3261SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr); 5869f074e33SMatthew G Knepley PetscFunctionReturn(0); 5879f074e33SMatthew G Knepley } 58807b0a1fcSMatthew G Knepley 58980330477SMatthew G. Knepley /*@ 59080330477SMatthew G. Knepley DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices 59180330477SMatthew G. Knepley 59280330477SMatthew G. Knepley Collective on DM 59380330477SMatthew G. Knepley 59480330477SMatthew G. Knepley Input Parameter: 59580330477SMatthew G. Knepley . dmA - The DMPlex object with initial coordinates 59680330477SMatthew G. Knepley 59780330477SMatthew G. Knepley Output Parameter: 59880330477SMatthew G. Knepley . dmB - The DMPlex object with copied coordinates 59980330477SMatthew G. Knepley 60080330477SMatthew G. Knepley Level: intermediate 60180330477SMatthew G. Knepley 60280330477SMatthew G. Knepley Note: This is typically used when adding pieces other than vertices to a mesh 60380330477SMatthew G. Knepley 60480330477SMatthew G. Knepley .keywords: mesh 60565f90189SJed Brown .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMGetCoordinateSection() 60680330477SMatthew G. Knepley @*/ 60707b0a1fcSMatthew G Knepley PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB) 60807b0a1fcSMatthew G Knepley { 60907b0a1fcSMatthew G Knepley Vec coordinatesA, coordinatesB; 61043eeeb2dSStefano Zampini VecType vtype; 61107b0a1fcSMatthew G Knepley PetscSection coordSectionA, coordSectionB; 61207b0a1fcSMatthew G Knepley PetscScalar *coordsA, *coordsB; 6130bedd6beSMatthew G. Knepley PetscInt spaceDim, Nf, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d; 61443eeeb2dSStefano Zampini PetscInt cStartA, cEndA, cStartB, cEndB, cS, cE; 61543eeeb2dSStefano Zampini PetscBool lc = PETSC_FALSE; 61607b0a1fcSMatthew G Knepley PetscErrorCode ierr; 61707b0a1fcSMatthew G Knepley 61807b0a1fcSMatthew G Knepley PetscFunctionBegin; 61943eeeb2dSStefano Zampini PetscValidHeaderSpecific(dmA, DM_CLASSID, 1); 62043eeeb2dSStefano Zampini PetscValidHeaderSpecific(dmB, DM_CLASSID, 2); 62176b791aaSMatthew G. Knepley if (dmA == dmB) PetscFunctionReturn(0); 62207b0a1fcSMatthew G Knepley ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr); 62307b0a1fcSMatthew G Knepley ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr); 62407b0a1fcSMatthew G Knepley if ((vEndA-vStartA) != (vEndB-vStartB)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of vertices in first DM %d != %d in the second DM", vEndA-vStartA, vEndB-vStartB); 62543eeeb2dSStefano Zampini ierr = DMPlexGetHeightStratum(dmA, 0, &cStartA, &cEndA);CHKERRQ(ierr); 62643eeeb2dSStefano Zampini ierr = DMPlexGetHeightStratum(dmB, 0, &cStartB, &cEndB);CHKERRQ(ierr); 62769d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr); 62869d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr); 629972bc18aSToby Isaac if (coordSectionA == coordSectionB) PetscFunctionReturn(0); 6300bedd6beSMatthew G. Knepley ierr = PetscSectionGetNumFields(coordSectionA, &Nf);CHKERRQ(ierr); 6310bedd6beSMatthew G. Knepley if (!Nf) PetscFunctionReturn(0); 6320bedd6beSMatthew G. Knepley if (Nf > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of coordinate fields must be 1, not %D", Nf); 633df26b574SMatthew G. Knepley if (!coordSectionB) { 634df26b574SMatthew G. Knepley PetscInt dim; 635df26b574SMatthew G. Knepley 636df26b574SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) coordSectionA), &coordSectionB);CHKERRQ(ierr); 637df26b574SMatthew G. Knepley ierr = DMGetCoordinateDim(dmA, &dim);CHKERRQ(ierr); 638df26b574SMatthew G. Knepley ierr = DMSetCoordinateSection(dmB, dim, coordSectionB);CHKERRQ(ierr); 639df26b574SMatthew G. Knepley ierr = PetscObjectDereference((PetscObject) coordSectionB);CHKERRQ(ierr); 640df26b574SMatthew G. Knepley } 64107b0a1fcSMatthew G Knepley ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr); 64207b0a1fcSMatthew G Knepley ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr); 64307b0a1fcSMatthew G Knepley ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr); 64443eeeb2dSStefano Zampini ierr = PetscSectionGetChart(coordSectionA, &cS, &cE);CHKERRQ(ierr); 64543eeeb2dSStefano Zampini if (cStartA <= cS && cS < cEndA) { /* localized coordinates */ 64643eeeb2dSStefano Zampini if ((cEndA-cStartA) != (cEndB-cStartB)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of cellls in first DM %d != %d in the second DM", cEndA-cStartA, cEndB-cStartB); 64743eeeb2dSStefano Zampini cS = cS - cStartA + cStartB; 64843eeeb2dSStefano Zampini cE = vEndB; 64943eeeb2dSStefano Zampini lc = PETSC_TRUE; 65043eeeb2dSStefano Zampini } else { 65143eeeb2dSStefano Zampini cS = vStartB; 65243eeeb2dSStefano Zampini cE = vEndB; 65343eeeb2dSStefano Zampini } 65443eeeb2dSStefano Zampini ierr = PetscSectionSetChart(coordSectionB, cS, cE);CHKERRQ(ierr); 65507b0a1fcSMatthew G Knepley for (v = vStartB; v < vEndB; ++v) { 65607b0a1fcSMatthew G Knepley ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr); 65707b0a1fcSMatthew G Knepley ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr); 65807b0a1fcSMatthew G Knepley } 65943eeeb2dSStefano Zampini if (lc) { /* localized coordinates */ 66043eeeb2dSStefano Zampini PetscInt c; 66143eeeb2dSStefano Zampini 66243eeeb2dSStefano Zampini for (c = cS-cStartB; c < cEndB-cStartB; c++) { 66343eeeb2dSStefano Zampini PetscInt dof; 66443eeeb2dSStefano Zampini 66543eeeb2dSStefano Zampini ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr); 66643eeeb2dSStefano Zampini ierr = PetscSectionSetDof(coordSectionB, c + cStartB, dof);CHKERRQ(ierr); 66743eeeb2dSStefano Zampini ierr = PetscSectionSetFieldDof(coordSectionB, c + cStartB, 0, dof);CHKERRQ(ierr); 66843eeeb2dSStefano Zampini } 66943eeeb2dSStefano Zampini } 67007b0a1fcSMatthew G Knepley ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr); 67107b0a1fcSMatthew G Knepley ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr); 67207b0a1fcSMatthew G Knepley ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr); 6738b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinatesB);CHKERRQ(ierr); 67407b0a1fcSMatthew G Knepley ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr); 67507b0a1fcSMatthew G Knepley ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr); 67643eeeb2dSStefano Zampini ierr = VecGetBlockSize(coordinatesA, &d);CHKERRQ(ierr); 67743eeeb2dSStefano Zampini ierr = VecSetBlockSize(coordinatesB, d);CHKERRQ(ierr); 67843eeeb2dSStefano Zampini ierr = VecGetType(coordinatesA, &vtype);CHKERRQ(ierr); 67943eeeb2dSStefano Zampini ierr = VecSetType(coordinatesB, vtype);CHKERRQ(ierr); 68007b0a1fcSMatthew G Knepley ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr); 68107b0a1fcSMatthew G Knepley ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr); 68207b0a1fcSMatthew G Knepley for (v = 0; v < vEndB-vStartB; ++v) { 68343eeeb2dSStefano Zampini PetscInt offA, offB; 68443eeeb2dSStefano Zampini 68543eeeb2dSStefano Zampini ierr = PetscSectionGetOffset(coordSectionA, v + vStartA, &offA);CHKERRQ(ierr); 68643eeeb2dSStefano Zampini ierr = PetscSectionGetOffset(coordSectionB, v + vStartB, &offB);CHKERRQ(ierr); 68707b0a1fcSMatthew G Knepley for (d = 0; d < spaceDim; ++d) { 68843eeeb2dSStefano Zampini coordsB[offB+d] = coordsA[offA+d]; 68943eeeb2dSStefano Zampini } 69043eeeb2dSStefano Zampini } 69143eeeb2dSStefano Zampini if (lc) { /* localized coordinates */ 69243eeeb2dSStefano Zampini PetscInt c; 69343eeeb2dSStefano Zampini 69443eeeb2dSStefano Zampini for (c = cS-cStartB; c < cEndB-cStartB; c++) { 69543eeeb2dSStefano Zampini PetscInt dof, offA, offB; 69643eeeb2dSStefano Zampini 69743eeeb2dSStefano Zampini ierr = PetscSectionGetOffset(coordSectionA, c + cStartA, &offA);CHKERRQ(ierr); 69843eeeb2dSStefano Zampini ierr = PetscSectionGetOffset(coordSectionB, c + cStartB, &offB);CHKERRQ(ierr); 69943eeeb2dSStefano Zampini ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr); 70043eeeb2dSStefano Zampini ierr = PetscMemcpy(coordsB + offB,coordsA + offA,dof*sizeof(*coordsB));CHKERRQ(ierr); 70107b0a1fcSMatthew G Knepley } 70207b0a1fcSMatthew G Knepley } 70307b0a1fcSMatthew G Knepley ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr); 70407b0a1fcSMatthew G Knepley ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr); 70507b0a1fcSMatthew G Knepley ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr); 70607b0a1fcSMatthew G Knepley ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr); 70707b0a1fcSMatthew G Knepley PetscFunctionReturn(0); 70807b0a1fcSMatthew G Knepley } 7095c386225SMatthew G. Knepley 7104e3744c5SMatthew G. Knepley /*@ 7114e3744c5SMatthew G. Knepley DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh 7124e3744c5SMatthew G. Knepley 7134e3744c5SMatthew G. Knepley Collective on DM 7144e3744c5SMatthew G. Knepley 7154e3744c5SMatthew G. Knepley Input Parameter: 7164e3744c5SMatthew G. Knepley . dm - The complete DMPlex object 7174e3744c5SMatthew G. Knepley 7184e3744c5SMatthew G. Knepley Output Parameter: 7194e3744c5SMatthew G. Knepley . dmUnint - The DMPlex object with only cells and vertices 7204e3744c5SMatthew G. Knepley 7214e3744c5SMatthew G. Knepley Level: intermediate 7224e3744c5SMatthew G. Knepley 72395452b02SPatrick Sanan Notes: 72495452b02SPatrick Sanan It does not copy over the coordinates. 72543eeeb2dSStefano Zampini 7264e3744c5SMatthew G. Knepley .keywords: mesh 72743eeeb2dSStefano Zampini .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates() 7284e3744c5SMatthew G. Knepley @*/ 7294e3744c5SMatthew G. Knepley PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint) 7304e3744c5SMatthew G. Knepley { 7314e3744c5SMatthew G. Knepley DM udm; 732595d4abbSMatthew G. Knepley PetscInt dim, vStart, vEnd, cStart, cEnd, c, maxConeSize = 0, *cone; 7334e3744c5SMatthew G. Knepley PetscErrorCode ierr; 7344e3744c5SMatthew G. Knepley 7354e3744c5SMatthew G. Knepley PetscFunctionBegin; 73643eeeb2dSStefano Zampini PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 73743eeeb2dSStefano Zampini PetscValidPointer(dmUnint, 2); 738c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 7394e3744c5SMatthew G. Knepley if (dim <= 1) { 7404e3744c5SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 741595d4abbSMatthew G. Knepley *dmUnint = dm; 742595d4abbSMatthew G. Knepley PetscFunctionReturn(0); 7434e3744c5SMatthew G. Knepley } 744595d4abbSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 745595d4abbSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 7464e3744c5SMatthew G. Knepley ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr); 7474e3744c5SMatthew G. Knepley ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr); 748c73cfb54SMatthew G. Knepley ierr = DMSetDimension(udm, dim);CHKERRQ(ierr); 7494e3744c5SMatthew G. Knepley ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr); 7504e3744c5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 751595d4abbSMatthew G. Knepley PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 7524e3744c5SMatthew G. Knepley 7534e3744c5SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 7544e3744c5SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 7554e3744c5SMatthew G. Knepley const PetscInt p = closure[cl]; 7564e3744c5SMatthew G. Knepley 7574e3744c5SMatthew G. Knepley if ((p >= vStart) && (p < vEnd)) ++coneSize; 7584e3744c5SMatthew G. Knepley } 7594e3744c5SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 7604e3744c5SMatthew G. Knepley ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr); 761595d4abbSMatthew G. Knepley maxConeSize = PetscMax(maxConeSize, coneSize); 7624e3744c5SMatthew G. Knepley } 7634e3744c5SMatthew G. Knepley ierr = DMSetUp(udm);CHKERRQ(ierr); 764785e854fSJed Brown ierr = PetscMalloc1(maxConeSize, &cone);CHKERRQ(ierr); 7654e3744c5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 766595d4abbSMatthew G. Knepley PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 7674e3744c5SMatthew G. Knepley 7684e3744c5SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 7694e3744c5SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 7704e3744c5SMatthew G. Knepley const PetscInt p = closure[cl]; 7714e3744c5SMatthew G. Knepley 7724e3744c5SMatthew G. Knepley if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p; 7734e3744c5SMatthew G. Knepley } 7744e3744c5SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 7754e3744c5SMatthew G. Knepley ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr); 7764e3744c5SMatthew G. Knepley } 7774e3744c5SMatthew G. Knepley ierr = PetscFree(cone);CHKERRQ(ierr); 7784e3744c5SMatthew G. Knepley ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr); 7794e3744c5SMatthew G. Knepley ierr = DMPlexStratify(udm);CHKERRQ(ierr); 7805c7de58aSMatthew G. Knepley /* Reduce SF */ 7815c7de58aSMatthew G. Knepley { 7825c7de58aSMatthew G. Knepley PetscSF sfPoint, sfPointUn; 7835c7de58aSMatthew G. Knepley const PetscSFNode *remotePoints; 7845c7de58aSMatthew G. Knepley const PetscInt *localPoints; 7855c7de58aSMatthew G. Knepley PetscSFNode *remotePointsUn; 7865c7de58aSMatthew G. Knepley PetscInt *localPointsUn; 7875c7de58aSMatthew G. Knepley PetscInt vEnd, numRoots, numLeaves, l; 7885c7de58aSMatthew G. Knepley PetscInt numLeavesUn = 0, n = 0; 7895c7de58aSMatthew G. Knepley PetscErrorCode ierr; 7905c7de58aSMatthew G. Knepley 7915c7de58aSMatthew G. Knepley /* Get original SF information */ 7925c7de58aSMatthew G. Knepley ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 7935c7de58aSMatthew G. Knepley ierr = DMGetPointSF(udm, &sfPointUn);CHKERRQ(ierr); 7945c7de58aSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, NULL, &vEnd);CHKERRQ(ierr); 7955c7de58aSMatthew G. Knepley ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 7965c7de58aSMatthew G. Knepley /* Allocate space for cells and vertices */ 7975c7de58aSMatthew G. Knepley for (l = 0; l < numLeaves; ++l) if (localPoints[l] < vEnd) numLeavesUn++; 7985c7de58aSMatthew G. Knepley /* Fill in leaves */ 7995c7de58aSMatthew G. Knepley if (vEnd >= 0) { 8005c7de58aSMatthew G. Knepley ierr = PetscMalloc1(numLeavesUn, &remotePointsUn);CHKERRQ(ierr); 8015c7de58aSMatthew G. Knepley ierr = PetscMalloc1(numLeavesUn, &localPointsUn);CHKERRQ(ierr); 8025c7de58aSMatthew G. Knepley for (l = 0; l < numLeaves; l++) { 8035c7de58aSMatthew G. Knepley if (localPoints[l] < vEnd) { 8045c7de58aSMatthew G. Knepley localPointsUn[n] = localPoints[l]; 8055c7de58aSMatthew G. Knepley remotePointsUn[n].rank = remotePoints[l].rank; 8065c7de58aSMatthew G. Knepley remotePointsUn[n].index = remotePoints[l].index; 8075c7de58aSMatthew G. Knepley ++n; 8085c7de58aSMatthew G. Knepley } 8095c7de58aSMatthew G. Knepley } 8105c7de58aSMatthew G. Knepley if (n != numLeavesUn) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent number of leaves %d != %d", n, numLeavesUn); 8115c7de58aSMatthew G. Knepley ierr = PetscSFSetGraph(sfPointUn, vEnd, numLeavesUn, localPointsUn, PETSC_OWN_POINTER, remotePointsUn, PETSC_OWN_POINTER);CHKERRQ(ierr); 8125c7de58aSMatthew G. Knepley } 8135c7de58aSMatthew G. Knepley } 81443eeeb2dSStefano Zampini { 81543eeeb2dSStefano Zampini PetscBool isper; 81643eeeb2dSStefano Zampini const PetscReal *maxCell, *L; 81743eeeb2dSStefano Zampini const DMBoundaryType *bd; 81843eeeb2dSStefano Zampini 81943eeeb2dSStefano Zampini ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr); 82043eeeb2dSStefano Zampini ierr = DMSetPeriodicity(udm,isper,maxCell,L,bd);CHKERRQ(ierr); 82143eeeb2dSStefano Zampini } 82243eeeb2dSStefano Zampini 8234e3744c5SMatthew G. Knepley *dmUnint = udm; 8244e3744c5SMatthew G. Knepley PetscFunctionReturn(0); 8254e3744c5SMatthew G. Knepley } 826