19f074e33SMatthew G Knepley #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 29f074e33SMatthew G Knepley #include <../src/sys/utils/hash.h> 39f074e33SMatthew G Knepley 49f074e33SMatthew G Knepley #undef __FUNCT__ 59a5b13a2SMatthew G Knepley #define __FUNCT__ "DMPlexGetFaces_Internal" 69f074e33SMatthew G Knepley /* 79a5b13a2SMatthew G Knepley DMPlexGetFaces_Internal - Gets groups of vertices that correspond to faces for the given cell 89f074e33SMatthew G Knepley */ 9439ece16SMatthew G. Knepley PetscErrorCode DMPlexGetFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 109f074e33SMatthew G Knepley { 119f074e33SMatthew G Knepley const PetscInt *cone = NULL; 129a5b13a2SMatthew G Knepley PetscInt maxConeSize, maxSupportSize, coneSize; 139f074e33SMatthew G Knepley PetscErrorCode ierr; 149f074e33SMatthew G Knepley 159f074e33SMatthew G Knepley PetscFunctionBegin; 169f074e33SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 179a5b13a2SMatthew G Knepley ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 189f074e33SMatthew G Knepley ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 199f074e33SMatthew G Knepley ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 20439ece16SMatthew G. Knepley ierr = DMPlexGetRawFaces_Internal(dm, dim, coneSize, cone, numFaces, faceSize, faces);CHKERRQ(ierr); 21439ece16SMatthew G. Knepley PetscFunctionReturn(0); 22439ece16SMatthew G. Knepley } 23439ece16SMatthew G. Knepley 24439ece16SMatthew G. Knepley #undef __FUNCT__ 25439ece16SMatthew G. Knepley #define __FUNCT__ "DMPlexRestoreFaces_Internal" 26439ece16SMatthew G. Knepley /* 27439ece16SMatthew G. Knepley DMPlexRestoreFaces_Internal - Restores the array 28439ece16SMatthew G. Knepley */ 29439ece16SMatthew G. Knepley PetscErrorCode DMPlexRestoreFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 30439ece16SMatthew G. Knepley { 31439ece16SMatthew G. Knepley PetscErrorCode ierr; 32439ece16SMatthew G. Knepley 33439ece16SMatthew G. Knepley PetscFunctionBegin; 3406acb219SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, 0, PETSC_INT, (void *) faces);CHKERRQ(ierr); 35439ece16SMatthew G. Knepley PetscFunctionReturn(0); 36439ece16SMatthew G. Knepley } 37439ece16SMatthew G. Knepley 38439ece16SMatthew G. Knepley #undef __FUNCT__ 39439ece16SMatthew G. Knepley #define __FUNCT__ "DMPlexGetRawFaces_Internal" 40439ece16SMatthew G. Knepley /* 41439ece16SMatthew G. Knepley DMPlexGetRawFaces_Internal - Gets groups of vertices that correspond to faces for the given cone 42439ece16SMatthew G. Knepley */ 43439ece16SMatthew G. Knepley PetscErrorCode DMPlexGetRawFaces_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 44439ece16SMatthew G. Knepley { 45439ece16SMatthew G. Knepley PetscInt *facesTmp; 46439ece16SMatthew G. Knepley PetscInt maxConeSize, maxSupportSize; 47439ece16SMatthew G. Knepley PetscErrorCode ierr; 48439ece16SMatthew G. Knepley 49439ece16SMatthew G. Knepley PetscFunctionBegin; 50439ece16SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 51439ece16SMatthew G. Knepley ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 528496a43fSMatthew G. Knepley if (faces) {ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), PETSC_INT, &facesTmp);CHKERRQ(ierr);} 539f074e33SMatthew G Knepley switch (dim) { 54*c49d9212SMatthew G. Knepley case 1: 55*c49d9212SMatthew G. Knepley switch (coneSize) { 56*c49d9212SMatthew G. Knepley case 2: 57*c49d9212SMatthew G. Knepley if (faces) { 58*c49d9212SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 59*c49d9212SMatthew G. Knepley *faces = facesTmp; 60*c49d9212SMatthew G. Knepley } 61*c49d9212SMatthew G. Knepley if (numFaces) *numFaces = 2; 62*c49d9212SMatthew G. Knepley if (faceSize) *faceSize = 1; 63*c49d9212SMatthew G. Knepley break; 64*c49d9212SMatthew G. Knepley default: 65*c49d9212SMatthew G. Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 66*c49d9212SMatthew G. Knepley } 67*c49d9212SMatthew G. Knepley break; 689f074e33SMatthew G Knepley case 2: 699f074e33SMatthew G Knepley switch (coneSize) { 709f074e33SMatthew G Knepley case 3: 719a5b13a2SMatthew G Knepley if (faces) { 729a5b13a2SMatthew G Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 739a5b13a2SMatthew G Knepley facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 749a5b13a2SMatthew G Knepley facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 759a5b13a2SMatthew G Knepley *faces = facesTmp; 769a5b13a2SMatthew G Knepley } 779a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 3; 789a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 2; 799f074e33SMatthew G Knepley break; 809f074e33SMatthew G Knepley case 4: 819a5b13a2SMatthew G Knepley /* Vertices follow right hand rule */ 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[3]; 869a5b13a2SMatthew G Knepley facesTmp[6] = cone[3]; facesTmp[7] = cone[0]; 879a5b13a2SMatthew G Knepley *faces = facesTmp; 889a5b13a2SMatthew G Knepley } 899a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 4; 909a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 2; 919f074e33SMatthew G Knepley break; 929f074e33SMatthew G Knepley default: 939f074e33SMatthew G Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 949f074e33SMatthew G Knepley } 959f074e33SMatthew G Knepley break; 969f074e33SMatthew G Knepley case 3: 979f074e33SMatthew G Knepley switch (coneSize) { 989f074e33SMatthew G Knepley case 3: 999a5b13a2SMatthew G Knepley if (faces) { 1009a5b13a2SMatthew G Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 1019a5b13a2SMatthew G Knepley facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 1029a5b13a2SMatthew G Knepley facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 1039a5b13a2SMatthew G Knepley *faces = facesTmp; 1049a5b13a2SMatthew G Knepley } 1059a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 3; 1069a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 2; 1079f074e33SMatthew G Knepley break; 1089f074e33SMatthew G Knepley case 4: 1092e1b13c2SMatthew G. Knepley /* Vertices of first face follow right hand rule and normal points away from last vertex */ 1109a5b13a2SMatthew G Knepley if (faces) { 1112e1b13c2SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; 1122e1b13c2SMatthew G. Knepley facesTmp[3] = cone[0]; facesTmp[4] = cone[3]; facesTmp[5] = cone[1]; 1132e1b13c2SMatthew G. Knepley facesTmp[6] = cone[0]; facesTmp[7] = cone[2]; facesTmp[8] = cone[3]; 1142e1b13c2SMatthew G. Knepley facesTmp[9] = cone[2]; facesTmp[10] = cone[1]; facesTmp[11] = cone[3]; 1159a5b13a2SMatthew G Knepley *faces = facesTmp; 1169a5b13a2SMatthew G Knepley } 1179a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 4; 1189a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 3; 1199a5b13a2SMatthew G Knepley break; 1209a5b13a2SMatthew G Knepley case 8: 1219a5b13a2SMatthew G Knepley if (faces) { 12251cfd6a4SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; facesTmp[3] = cone[3]; /* Bottom */ 12351cfd6a4SMatthew G. Knepley facesTmp[4] = cone[4]; facesTmp[5] = cone[5]; facesTmp[6] = cone[6]; facesTmp[7] = cone[7]; /* Top */ 12451cfd6a4SMatthew G. Knepley facesTmp[8] = cone[0]; facesTmp[9] = cone[3]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4]; /* Front */ 12551cfd6a4SMatthew G. Knepley facesTmp[12] = cone[2]; facesTmp[13] = cone[1]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6]; /* Back */ 12651cfd6a4SMatthew G. Knepley facesTmp[16] = cone[3]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5]; /* Right */ 12751cfd6a4SMatthew G. Knepley facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[1]; /* Left */ 1289a5b13a2SMatthew G Knepley *faces = facesTmp; 1299a5b13a2SMatthew G Knepley } 1309a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 6; 1319a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 4; 1329f074e33SMatthew G Knepley break; 1339f074e33SMatthew G Knepley default: 1349f074e33SMatthew G Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 1359f074e33SMatthew G Knepley } 1369f074e33SMatthew G Knepley break; 1379f074e33SMatthew G Knepley default: 1389f074e33SMatthew G Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim); 1399f074e33SMatthew G Knepley } 1409f074e33SMatthew G Knepley PetscFunctionReturn(0); 1419f074e33SMatthew G Knepley } 1429f074e33SMatthew G Knepley 1439f074e33SMatthew G Knepley #undef __FUNCT__ 1449a5b13a2SMatthew G Knepley #define __FUNCT__ "DMPlexInterpolateFaces_Internal" 1459a5b13a2SMatthew G Knepley /* This interpolates faces for cells at some stratum */ 1469a5b13a2SMatthew G Knepley static PetscErrorCode DMPlexInterpolateFaces_Internal(DM dm, PetscInt cellDepth, DM idm) 1479f074e33SMatthew G Knepley { 148d4efc6f1SMatthew G. Knepley DMLabel subpointMap; 1499a5b13a2SMatthew G Knepley PetscHashIJKL faceTable; 1509a5b13a2SMatthew G Knepley PetscInt *pStart, *pEnd; 1519a5b13a2SMatthew G Knepley PetscInt cellDim, depth, faceDepth = cellDepth, numPoints = 0, faceSizeAll = 0, face, c, d; 1529f074e33SMatthew G Knepley PetscErrorCode ierr; 1539f074e33SMatthew G Knepley 1549f074e33SMatthew G Knepley PetscFunctionBegin; 1559a5b13a2SMatthew G Knepley ierr = DMPlexGetDimension(dm, &cellDim);CHKERRQ(ierr); 156d4efc6f1SMatthew G. Knepley /* HACK: I need a better way to determine face dimension, or an alternative to GetFaces() */ 157d4efc6f1SMatthew G. Knepley ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr); 158d4efc6f1SMatthew G. Knepley if (subpointMap) ++cellDim; 1599a5b13a2SMatthew G Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1609a5b13a2SMatthew G Knepley ++depth; 1619a5b13a2SMatthew G Knepley ++cellDepth; 1629a5b13a2SMatthew G Knepley cellDim -= depth - cellDepth; 1639a5b13a2SMatthew G Knepley ierr = PetscMalloc2(depth+1,PetscInt,&pStart,depth+1,PetscInt,&pEnd);CHKERRQ(ierr); 1649a5b13a2SMatthew G Knepley for (d = depth-1; d >= faceDepth; --d) { 1659a5b13a2SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, d, &pStart[d+1], &pEnd[d+1]);CHKERRQ(ierr); 1669f074e33SMatthew G Knepley } 1679a5b13a2SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, -1, NULL, &pStart[faceDepth]);CHKERRQ(ierr); 1689a5b13a2SMatthew G Knepley pEnd[faceDepth] = pStart[faceDepth]; 1699a5b13a2SMatthew G Knepley for (d = faceDepth-1; d >= 0; --d) { 1709a5b13a2SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 1719f074e33SMatthew G Knepley } 1729a5b13a2SMatthew G Knepley if (pEnd[cellDepth] > pStart[cellDepth]) {ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr);} 1739a5b13a2SMatthew 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); 1749a5b13a2SMatthew G Knepley ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 1759a5b13a2SMatthew G Knepley ierr = PetscHashIJKLSetMultivalued(faceTable, PETSC_FALSE);CHKERRQ(ierr); 1769a5b13a2SMatthew G Knepley for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) { 1779f074e33SMatthew G Knepley const PetscInt *cellFaces; 1789a5b13a2SMatthew G Knepley PetscInt numCellFaces, faceSize, cf, f; 1799f074e33SMatthew G Knepley 1809a5b13a2SMatthew G Knepley ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 1819a5b13a2SMatthew G Knepley if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 1829f074e33SMatthew G Knepley for (cf = 0; cf < numCellFaces; ++cf) { 1839a5b13a2SMatthew G Knepley const PetscInt *cellFace = &cellFaces[cf*faceSize]; 1849a5b13a2SMatthew G Knepley PetscHashIJKLKey key; 1859f074e33SMatthew G Knepley 1869a5b13a2SMatthew G Knepley if (faceSize == 2) { 1879a5b13a2SMatthew G Knepley key.i = PetscMin(cellFace[0], cellFace[1]); 1889a5b13a2SMatthew G Knepley key.j = PetscMax(cellFace[0], cellFace[1]); 1897fcd4ef0SJed Brown key.k = 0; 1907fcd4ef0SJed Brown key.l = 0; 1919a5b13a2SMatthew G Knepley } else { 1929a5b13a2SMatthew G Knepley key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0; 1939a5b13a2SMatthew G Knepley ierr = PetscSortInt(faceSize, (PetscInt *) &key); 1949f074e33SMatthew G Knepley } 1959a5b13a2SMatthew G Knepley ierr = PetscHashIJKLGet(faceTable, key, &f);CHKERRQ(ierr); 1969a5b13a2SMatthew G Knepley if (f < 0) { 1979a5b13a2SMatthew G Knepley ierr = PetscHashIJKLAdd(faceTable, key, face);CHKERRQ(ierr); 1989a5b13a2SMatthew G Knepley f = face++; 1999a5b13a2SMatthew G Knepley } 2009a5b13a2SMatthew G Knepley } 201439ece16SMatthew G. Knepley ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 2029a5b13a2SMatthew G Knepley } 2039a5b13a2SMatthew G Knepley pEnd[faceDepth] = face; 2049a5b13a2SMatthew G Knepley ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 2059a5b13a2SMatthew G Knepley /* Count new points */ 2069a5b13a2SMatthew G Knepley for (d = 0; d <= depth; ++d) { 2079a5b13a2SMatthew G Knepley numPoints += pEnd[d]-pStart[d]; 2089a5b13a2SMatthew G Knepley } 2099a5b13a2SMatthew G Knepley ierr = DMPlexSetChart(idm, 0, numPoints);CHKERRQ(ierr); 2109a5b13a2SMatthew G Knepley /* Set cone sizes */ 2119a5b13a2SMatthew G Knepley for (d = 0; d <= depth; ++d) { 2129a5b13a2SMatthew G Knepley PetscInt coneSize, p; 2139f074e33SMatthew G Knepley 2149a5b13a2SMatthew G Knepley if (d == faceDepth) { 2159a5b13a2SMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 2169a5b13a2SMatthew G Knepley /* I see no way to do this if we admit faces of different shapes */ 2179a5b13a2SMatthew G Knepley ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr); 2189a5b13a2SMatthew G Knepley } 219a014044eSMatthew G Knepley } else if (d == cellDepth) { 220a014044eSMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 221a014044eSMatthew G Knepley /* Number of cell faces may be different from number of cell vertices*/ 222a014044eSMatthew G Knepley ierr = DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);CHKERRQ(ierr); 223a014044eSMatthew G Knepley ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 224a014044eSMatthew G Knepley } 2259a5b13a2SMatthew G Knepley } else { 2269a5b13a2SMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 2279a5b13a2SMatthew G Knepley ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 2289a5b13a2SMatthew G Knepley ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 2299f074e33SMatthew G Knepley } 2309f074e33SMatthew G Knepley } 2319f074e33SMatthew G Knepley } 2329f074e33SMatthew G Knepley ierr = DMSetUp(idm);CHKERRQ(ierr); 2339f074e33SMatthew G Knepley /* Get face cones from subsets of cell vertices */ 2349a5b13a2SMatthew 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); 2359a5b13a2SMatthew G Knepley ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 2369a5b13a2SMatthew G Knepley ierr = PetscHashIJKLSetMultivalued(faceTable, PETSC_FALSE);CHKERRQ(ierr); 2379a5b13a2SMatthew G Knepley for (d = depth; d > cellDepth; --d) { 2389a5b13a2SMatthew G Knepley const PetscInt *cone; 2399a5b13a2SMatthew G Knepley PetscInt p; 2409a5b13a2SMatthew G Knepley 2419a5b13a2SMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 2429a5b13a2SMatthew G Knepley ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 2439a5b13a2SMatthew G Knepley ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr); 2449a5b13a2SMatthew G Knepley ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr); 2459a5b13a2SMatthew G Knepley ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr); 2469f074e33SMatthew G Knepley } 2479a5b13a2SMatthew G Knepley } 2489a5b13a2SMatthew G Knepley for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) { 2499f074e33SMatthew G Knepley const PetscInt *cellFaces; 2509a5b13a2SMatthew G Knepley PetscInt numCellFaces, faceSize, cf, f; 2519f074e33SMatthew G Knepley 2529a5b13a2SMatthew G Knepley ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 2539a5b13a2SMatthew G Knepley if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 2549f074e33SMatthew G Knepley for (cf = 0; cf < numCellFaces; ++cf) { 2559a5b13a2SMatthew G Knepley const PetscInt *cellFace = &cellFaces[cf*faceSize]; 2569a5b13a2SMatthew G Knepley PetscHashIJKLKey key; 2579f074e33SMatthew G Knepley 2589a5b13a2SMatthew G Knepley if (faceSize == 2) { 2599a5b13a2SMatthew G Knepley key.i = PetscMin(cellFace[0], cellFace[1]); 2609a5b13a2SMatthew G Knepley key.j = PetscMax(cellFace[0], cellFace[1]); 2617fcd4ef0SJed Brown key.k = 0; 2627fcd4ef0SJed Brown key.l = 0; 2639a5b13a2SMatthew G Knepley } else { 2649a5b13a2SMatthew G Knepley key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0; 2659a5b13a2SMatthew G Knepley ierr = PetscSortInt(faceSize, (PetscInt *) &key); 2669f074e33SMatthew G Knepley } 2679a5b13a2SMatthew G Knepley ierr = PetscHashIJKLGet(faceTable, key, &f);CHKERRQ(ierr); 2689a5b13a2SMatthew G Knepley if (f < 0) { 2699a5b13a2SMatthew G Knepley ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr); 2709a5b13a2SMatthew G Knepley ierr = PetscHashIJKLAdd(faceTable, key, face);CHKERRQ(ierr); 2719a5b13a2SMatthew G Knepley f = face++; 2729f074e33SMatthew G Knepley ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr); 2739a5b13a2SMatthew G Knepley } else { 2749a5b13a2SMatthew G Knepley const PetscInt *cone; 2759a5b13a2SMatthew G Knepley PetscInt coneSize, ornt, i, j; 2769f074e33SMatthew G Knepley 2779a5b13a2SMatthew G Knepley ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr); 2782e1b13c2SMatthew G. Knepley /* Orient face: Do not allow reverse orientation at the first vertex */ 2799f074e33SMatthew G Knepley ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr); 2809f074e33SMatthew G Knepley ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr); 2819a5b13a2SMatthew 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); 2829a5b13a2SMatthew G Knepley /* - First find the initial vertex */ 2839a5b13a2SMatthew G Knepley for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break; 2849a5b13a2SMatthew G Knepley /* - Try forward comparison */ 2859a5b13a2SMatthew G Knepley for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break; 2869a5b13a2SMatthew G Knepley if (j == faceSize) { 2879a5b13a2SMatthew G Knepley if ((faceSize == 2) && (i == 1)) ornt = -2; 2889a5b13a2SMatthew G Knepley else ornt = i; 2899a5b13a2SMatthew G Knepley } else { 2909a5b13a2SMatthew G Knepley /* - Try backward comparison */ 2919a5b13a2SMatthew G Knepley for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break; 2922e1b13c2SMatthew G. Knepley if (j == faceSize) { 2932e1b13c2SMatthew G. Knepley if (i == 0) ornt = -faceSize; 294dc1a705cSMatthew G. Knepley else ornt = -i; 2952e1b13c2SMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation"); 2969a5b13a2SMatthew G Knepley } 2979a5b13a2SMatthew G Knepley ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr); 2989f074e33SMatthew G Knepley } 2999f074e33SMatthew G Knepley } 300439ece16SMatthew G. Knepley ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 3019f074e33SMatthew G Knepley } 3029a5b13a2SMatthew 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]); 303c907b753SJed Brown ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 3049a5b13a2SMatthew G Knepley ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 3056551a8c7SMatthew G. Knepley ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 3069a5b13a2SMatthew G Knepley ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr); 3079a5b13a2SMatthew G Knepley ierr = DMPlexStratify(idm);CHKERRQ(ierr); 3089f074e33SMatthew G Knepley PetscFunctionReturn(0); 3099f074e33SMatthew G Knepley } 3109f074e33SMatthew G Knepley 3119f074e33SMatthew G Knepley #undef __FUNCT__ 3129f074e33SMatthew G Knepley #define __FUNCT__ "DMPlexInterpolate" 31380330477SMatthew G. Knepley /*@ 31480330477SMatthew G. Knepley DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc. 31580330477SMatthew G. Knepley 31680330477SMatthew G. Knepley Collective on DM 31780330477SMatthew G. Knepley 31880330477SMatthew G. Knepley Input Parameter: 3194e3744c5SMatthew G. Knepley . dm - The DMPlex object with only cells and vertices 32080330477SMatthew G. Knepley 32180330477SMatthew G. Knepley Output Parameter: 3224e3744c5SMatthew G. Knepley . dmInt - The complete DMPlex object 32380330477SMatthew G. Knepley 32480330477SMatthew G. Knepley Level: intermediate 32580330477SMatthew G. Knepley 32680330477SMatthew G. Knepley .keywords: mesh 3274e3744c5SMatthew G. Knepley .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList() 32880330477SMatthew G. Knepley @*/ 3299f074e33SMatthew G Knepley PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt) 3309f074e33SMatthew G Knepley { 3319a5b13a2SMatthew G Knepley DM idm, odm = dm; 3322e1b13c2SMatthew G. Knepley PetscInt depth, dim, d; 3339f074e33SMatthew G Knepley PetscErrorCode ierr; 3349f074e33SMatthew G Knepley 3359f074e33SMatthew G Knepley PetscFunctionBegin; 3362e1b13c2SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 3379f074e33SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 33876b791aaSMatthew G. Knepley if (dim <= 1) { 33976b791aaSMatthew G. Knepley ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 34076b791aaSMatthew G. Knepley idm = dm; 34176b791aaSMatthew G. Knepley } 3429a5b13a2SMatthew G Knepley for (d = 1; d < dim; ++d) { 3439a5b13a2SMatthew G Knepley /* Create interpolated mesh */ 3449a5b13a2SMatthew G Knepley ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr); 3459a5b13a2SMatthew G Knepley ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr); 3469a5b13a2SMatthew G Knepley ierr = DMPlexSetDimension(idm, dim);CHKERRQ(ierr); 3472e1b13c2SMatthew G. Knepley if (depth > 0) {ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr);} 3489a5b13a2SMatthew G Knepley if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);} 3499a5b13a2SMatthew G Knepley odm = idm; 3509f074e33SMatthew G Knepley } 3519a5b13a2SMatthew G Knepley *dmInt = idm; 3529f074e33SMatthew G Knepley PetscFunctionReturn(0); 3539f074e33SMatthew G Knepley } 35407b0a1fcSMatthew G Knepley 35507b0a1fcSMatthew G Knepley #undef __FUNCT__ 35607b0a1fcSMatthew G Knepley #define __FUNCT__ "DMPlexCopyCoordinates" 35780330477SMatthew G. Knepley /*@ 35880330477SMatthew G. Knepley DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices 35980330477SMatthew G. Knepley 36080330477SMatthew G. Knepley Collective on DM 36180330477SMatthew G. Knepley 36280330477SMatthew G. Knepley Input Parameter: 36380330477SMatthew G. Knepley . dmA - The DMPlex object with initial coordinates 36480330477SMatthew G. Knepley 36580330477SMatthew G. Knepley Output Parameter: 36680330477SMatthew G. Knepley . dmB - The DMPlex object with copied coordinates 36780330477SMatthew G. Knepley 36880330477SMatthew G. Knepley Level: intermediate 36980330477SMatthew G. Knepley 37080330477SMatthew G. Knepley Note: This is typically used when adding pieces other than vertices to a mesh 37180330477SMatthew G. Knepley 37280330477SMatthew G. Knepley .keywords: mesh 3735c386225SMatthew G. Knepley .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMPlexGetCoordinateSection() 37480330477SMatthew G. Knepley @*/ 37507b0a1fcSMatthew G Knepley PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB) 37607b0a1fcSMatthew G Knepley { 37707b0a1fcSMatthew G Knepley Vec coordinatesA, coordinatesB; 37807b0a1fcSMatthew G Knepley PetscSection coordSectionA, coordSectionB; 37907b0a1fcSMatthew G Knepley PetscScalar *coordsA, *coordsB; 38007b0a1fcSMatthew G Knepley PetscInt spaceDim, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d; 38107b0a1fcSMatthew G Knepley PetscErrorCode ierr; 38207b0a1fcSMatthew G Knepley 38307b0a1fcSMatthew G Knepley PetscFunctionBegin; 38476b791aaSMatthew G. Knepley if (dmA == dmB) PetscFunctionReturn(0); 38507b0a1fcSMatthew G Knepley ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr); 38607b0a1fcSMatthew G Knepley ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr); 38707b0a1fcSMatthew 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); 38807b0a1fcSMatthew G Knepley ierr = DMPlexGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr); 38907b0a1fcSMatthew G Knepley ierr = DMPlexGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr); 39007b0a1fcSMatthew G Knepley ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr); 39107b0a1fcSMatthew G Knepley ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr); 39207b0a1fcSMatthew G Knepley ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr); 39307b0a1fcSMatthew G Knepley ierr = PetscSectionSetChart(coordSectionB, vStartB, vEndB);CHKERRQ(ierr); 39407b0a1fcSMatthew G Knepley for (v = vStartB; v < vEndB; ++v) { 39507b0a1fcSMatthew G Knepley ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr); 39607b0a1fcSMatthew G Knepley ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr); 39707b0a1fcSMatthew G Knepley } 39807b0a1fcSMatthew G Knepley ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr); 39907b0a1fcSMatthew G Knepley ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr); 40007b0a1fcSMatthew G Knepley ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr); 40107b0a1fcSMatthew G Knepley ierr = VecCreate(PetscObjectComm((PetscObject) dmB), &coordinatesB);CHKERRQ(ierr); 40207b0a1fcSMatthew G Knepley ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr); 40307b0a1fcSMatthew G Knepley ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr); 4042eb5907fSJed Brown ierr = VecSetType(coordinatesB,VECSTANDARD);CHKERRQ(ierr); 40507b0a1fcSMatthew G Knepley ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr); 40607b0a1fcSMatthew G Knepley ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr); 40707b0a1fcSMatthew G Knepley for (v = 0; v < vEndB-vStartB; ++v) { 40807b0a1fcSMatthew G Knepley for (d = 0; d < spaceDim; ++d) { 40907b0a1fcSMatthew G Knepley coordsB[v*spaceDim+d] = coordsA[v*spaceDim+d]; 41007b0a1fcSMatthew G Knepley } 41107b0a1fcSMatthew G Knepley } 41207b0a1fcSMatthew G Knepley ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr); 41307b0a1fcSMatthew G Knepley ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr); 41407b0a1fcSMatthew G Knepley ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr); 41507b0a1fcSMatthew G Knepley ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr); 41607b0a1fcSMatthew G Knepley PetscFunctionReturn(0); 41707b0a1fcSMatthew G Knepley } 4185c386225SMatthew G. Knepley 4195c386225SMatthew G. Knepley #undef __FUNCT__ 4205c386225SMatthew G. Knepley #define __FUNCT__ "DMPlexCopyLabels" 4215c386225SMatthew G. Knepley /*@ 4225c386225SMatthew G. Knepley DMPlexCopyLabels - Copy labels from one mesh to another with a superset of the points 4235c386225SMatthew G. Knepley 4245c386225SMatthew G. Knepley Collective on DM 4255c386225SMatthew G. Knepley 4265c386225SMatthew G. Knepley Input Parameter: 4275c386225SMatthew G. Knepley . dmA - The DMPlex object with initial labels 4285c386225SMatthew G. Knepley 4295c386225SMatthew G. Knepley Output Parameter: 4305c386225SMatthew G. Knepley . dmB - The DMPlex object with copied labels 4315c386225SMatthew G. Knepley 4325c386225SMatthew G. Knepley Level: intermediate 4335c386225SMatthew G. Knepley 4345c386225SMatthew G. Knepley Note: This is typically used when interpolating or otherwise adding to a mesh 4355c386225SMatthew G. Knepley 4365c386225SMatthew G. Knepley .keywords: mesh 4375c386225SMatthew G. Knepley .seealso: DMCopyCoordinates(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMPlexGetCoordinateSection() 4385c386225SMatthew G. Knepley @*/ 4395c386225SMatthew G. Knepley PetscErrorCode DMPlexCopyLabels(DM dmA, DM dmB) 4405c386225SMatthew G. Knepley { 4415c386225SMatthew G. Knepley PetscInt numLabels, l; 4425c386225SMatthew G. Knepley PetscErrorCode ierr; 4435c386225SMatthew G. Knepley 4445c386225SMatthew G. Knepley PetscFunctionBegin; 4455c386225SMatthew G. Knepley if (dmA == dmB) PetscFunctionReturn(0); 4465c386225SMatthew G. Knepley ierr = DMPlexGetNumLabels(dmA, &numLabels);CHKERRQ(ierr); 4475c386225SMatthew G. Knepley for (l = 0; l < numLabels; ++l) { 4485c386225SMatthew G. Knepley DMLabel label, labelNew; 4495c386225SMatthew G. Knepley const char *name; 4505c386225SMatthew G. Knepley PetscBool flg; 4515c386225SMatthew G. Knepley 4525c386225SMatthew G. Knepley ierr = DMPlexGetLabelName(dmA, l, &name);CHKERRQ(ierr); 4535c386225SMatthew G. Knepley ierr = PetscStrcmp(name, "depth", &flg);CHKERRQ(ierr); 4545c386225SMatthew G. Knepley if (flg) continue; 4555c386225SMatthew G. Knepley ierr = DMPlexGetLabel(dmA, name, &label);CHKERRQ(ierr); 4565c386225SMatthew G. Knepley ierr = DMLabelDuplicate(label, &labelNew);CHKERRQ(ierr); 4575c386225SMatthew G. Knepley ierr = DMPlexAddLabel(dmB, labelNew);CHKERRQ(ierr); 4585c386225SMatthew G. Knepley } 4595c386225SMatthew G. Knepley PetscFunctionReturn(0); 4605c386225SMatthew G. Knepley } 461f5cedc29SMatthew G. Knepley 462f5cedc29SMatthew G. Knepley #undef __FUNCT__ 4634e3744c5SMatthew G. Knepley #define __FUNCT__ "DMPlexUninterpolate" 4644e3744c5SMatthew G. Knepley /*@ 4654e3744c5SMatthew G. Knepley DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh 4664e3744c5SMatthew G. Knepley 4674e3744c5SMatthew G. Knepley Collective on DM 4684e3744c5SMatthew G. Knepley 4694e3744c5SMatthew G. Knepley Input Parameter: 4704e3744c5SMatthew G. Knepley . dm - The complete DMPlex object 4714e3744c5SMatthew G. Knepley 4724e3744c5SMatthew G. Knepley Output Parameter: 4734e3744c5SMatthew G. Knepley . dmUnint - The DMPlex object with only cells and vertices 4744e3744c5SMatthew G. Knepley 4754e3744c5SMatthew G. Knepley Level: intermediate 4764e3744c5SMatthew G. Knepley 4774e3744c5SMatthew G. Knepley .keywords: mesh 4784e3744c5SMatthew G. Knepley .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList() 4794e3744c5SMatthew G. Knepley @*/ 4804e3744c5SMatthew G. Knepley PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint) 4814e3744c5SMatthew G. Knepley { 4824e3744c5SMatthew G. Knepley DM udm; 483595d4abbSMatthew G. Knepley PetscInt dim, vStart, vEnd, cStart, cEnd, c, maxConeSize = 0, *cone; 4844e3744c5SMatthew G. Knepley PetscErrorCode ierr; 4854e3744c5SMatthew G. Knepley 4864e3744c5SMatthew G. Knepley PetscFunctionBegin; 4874e3744c5SMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 4884e3744c5SMatthew G. Knepley if (dim <= 1) { 4894e3744c5SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 490595d4abbSMatthew G. Knepley *dmUnint = dm; 491595d4abbSMatthew G. Knepley PetscFunctionReturn(0); 4924e3744c5SMatthew G. Knepley } 493595d4abbSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 494595d4abbSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 4954e3744c5SMatthew G. Knepley ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr); 4964e3744c5SMatthew G. Knepley ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr); 4974e3744c5SMatthew G. Knepley ierr = DMPlexSetDimension(udm, dim);CHKERRQ(ierr); 4984e3744c5SMatthew G. Knepley ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr); 4994e3744c5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 500595d4abbSMatthew G. Knepley PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 5014e3744c5SMatthew G. Knepley 5024e3744c5SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 5034e3744c5SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 5044e3744c5SMatthew G. Knepley const PetscInt p = closure[cl]; 5054e3744c5SMatthew G. Knepley 5064e3744c5SMatthew G. Knepley if ((p >= vStart) && (p < vEnd)) ++coneSize; 5074e3744c5SMatthew G. Knepley } 5084e3744c5SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 5094e3744c5SMatthew G. Knepley ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr); 510595d4abbSMatthew G. Knepley maxConeSize = PetscMax(maxConeSize, coneSize); 5114e3744c5SMatthew G. Knepley } 5124e3744c5SMatthew G. Knepley ierr = DMSetUp(udm);CHKERRQ(ierr); 5134e3744c5SMatthew G. Knepley ierr = PetscMalloc(maxConeSize * sizeof(PetscInt), &cone);CHKERRQ(ierr); 5144e3744c5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 515595d4abbSMatthew G. Knepley PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 5164e3744c5SMatthew G. Knepley 5174e3744c5SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 5184e3744c5SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 5194e3744c5SMatthew G. Knepley const PetscInt p = closure[cl]; 5204e3744c5SMatthew G. Knepley 5214e3744c5SMatthew G. Knepley if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p; 5224e3744c5SMatthew G. Knepley } 5234e3744c5SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 5244e3744c5SMatthew G. Knepley ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr); 5254e3744c5SMatthew G. Knepley } 5264e3744c5SMatthew G. Knepley ierr = PetscFree(cone);CHKERRQ(ierr); 5274e3744c5SMatthew G. Knepley ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr); 5284e3744c5SMatthew G. Knepley ierr = DMPlexStratify(udm);CHKERRQ(ierr); 5294e3744c5SMatthew G. Knepley *dmUnint = udm; 5304e3744c5SMatthew G. Knepley PetscFunctionReturn(0); 5314e3744c5SMatthew G. Knepley } 532