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) { 549f074e33SMatthew G Knepley case 2: 559f074e33SMatthew G Knepley switch (coneSize) { 569f074e33SMatthew G Knepley case 3: 579a5b13a2SMatthew G Knepley if (faces) { 589a5b13a2SMatthew G Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 599a5b13a2SMatthew G Knepley facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 609a5b13a2SMatthew G Knepley facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 619a5b13a2SMatthew G Knepley *faces = facesTmp; 629a5b13a2SMatthew G Knepley } 639a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 3; 649a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 2; 659f074e33SMatthew G Knepley break; 669f074e33SMatthew G Knepley case 4: 679a5b13a2SMatthew G Knepley /* Vertices follow right hand rule */ 689a5b13a2SMatthew G Knepley if (faces) { 699a5b13a2SMatthew G Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 709a5b13a2SMatthew G Knepley facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 719a5b13a2SMatthew G Knepley facesTmp[4] = cone[2]; facesTmp[5] = cone[3]; 729a5b13a2SMatthew G Knepley facesTmp[6] = cone[3]; facesTmp[7] = cone[0]; 739a5b13a2SMatthew G Knepley *faces = facesTmp; 749a5b13a2SMatthew G Knepley } 759a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 4; 769a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 2; 779f074e33SMatthew G Knepley break; 789f074e33SMatthew G Knepley default: 799f074e33SMatthew G Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 809f074e33SMatthew G Knepley } 819f074e33SMatthew G Knepley break; 829f074e33SMatthew G Knepley case 3: 839f074e33SMatthew G Knepley switch (coneSize) { 849f074e33SMatthew G Knepley case 3: 859a5b13a2SMatthew G Knepley if (faces) { 869a5b13a2SMatthew G Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 879a5b13a2SMatthew G Knepley facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 889a5b13a2SMatthew G Knepley facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 899a5b13a2SMatthew G Knepley *faces = facesTmp; 909a5b13a2SMatthew G Knepley } 919a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 3; 929a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 2; 939f074e33SMatthew G Knepley break; 949f074e33SMatthew G Knepley case 4: 952e1b13c2SMatthew G. Knepley /* Vertices of first face follow right hand rule and normal points away from last vertex */ 969a5b13a2SMatthew G Knepley if (faces) { 972e1b13c2SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; 982e1b13c2SMatthew G. Knepley facesTmp[3] = cone[0]; facesTmp[4] = cone[3]; facesTmp[5] = cone[1]; 992e1b13c2SMatthew G. Knepley facesTmp[6] = cone[0]; facesTmp[7] = cone[2]; facesTmp[8] = cone[3]; 1002e1b13c2SMatthew G. Knepley facesTmp[9] = cone[2]; facesTmp[10] = cone[1]; facesTmp[11] = cone[3]; 1019a5b13a2SMatthew G Knepley *faces = facesTmp; 1029a5b13a2SMatthew G Knepley } 1039a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 4; 1049a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 3; 1059a5b13a2SMatthew G Knepley break; 1069a5b13a2SMatthew G Knepley case 8: 1079a5b13a2SMatthew G Knepley if (faces) { 1082e1b13c2SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; facesTmp[3] = cone[3]; 109a014044eSMatthew G Knepley facesTmp[4] = cone[4]; facesTmp[5] = cone[5]; facesTmp[6] = cone[6]; facesTmp[7] = cone[7]; 1102e1b13c2SMatthew G. Knepley facesTmp[8] = cone[0]; facesTmp[9] = cone[3]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4]; 1112e1b13c2SMatthew G. Knepley facesTmp[12] = cone[2]; facesTmp[13] = cone[1]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6]; 1122e1b13c2SMatthew G. Knepley facesTmp[16] = cone[3]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5]; 1132e1b13c2SMatthew G. Knepley facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[1]; 1149a5b13a2SMatthew G Knepley *faces = facesTmp; 1159a5b13a2SMatthew G Knepley } 1169a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 6; 1179a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 4; 1189f074e33SMatthew G Knepley break; 1199f074e33SMatthew G Knepley default: 1209f074e33SMatthew G Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 1219f074e33SMatthew G Knepley } 1229f074e33SMatthew G Knepley break; 1239f074e33SMatthew G Knepley default: 1249f074e33SMatthew G Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim); 1259f074e33SMatthew G Knepley } 1269f074e33SMatthew G Knepley PetscFunctionReturn(0); 1279f074e33SMatthew G Knepley } 1289f074e33SMatthew G Knepley 1299f074e33SMatthew G Knepley #undef __FUNCT__ 1309a5b13a2SMatthew G Knepley #define __FUNCT__ "DMPlexInterpolateFaces_Internal" 1319a5b13a2SMatthew G Knepley /* This interpolates faces for cells at some stratum */ 1329a5b13a2SMatthew G Knepley static PetscErrorCode DMPlexInterpolateFaces_Internal(DM dm, PetscInt cellDepth, DM idm) 1339f074e33SMatthew G Knepley { 134d4efc6f1SMatthew G. Knepley DMLabel subpointMap; 1359a5b13a2SMatthew G Knepley PetscHashIJKL faceTable; 1369a5b13a2SMatthew G Knepley PetscInt *pStart, *pEnd; 1379a5b13a2SMatthew G Knepley PetscInt cellDim, depth, faceDepth = cellDepth, numPoints = 0, faceSizeAll = 0, face, c, d; 1389f074e33SMatthew G Knepley PetscErrorCode ierr; 1399f074e33SMatthew G Knepley 1409f074e33SMatthew G Knepley PetscFunctionBegin; 1419a5b13a2SMatthew G Knepley ierr = DMPlexGetDimension(dm, &cellDim);CHKERRQ(ierr); 142d4efc6f1SMatthew G. Knepley /* HACK: I need a better way to determine face dimension, or an alternative to GetFaces() */ 143d4efc6f1SMatthew G. Knepley ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr); 144d4efc6f1SMatthew G. Knepley if (subpointMap) ++cellDim; 1459a5b13a2SMatthew G Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1469a5b13a2SMatthew G Knepley ++depth; 1479a5b13a2SMatthew G Knepley ++cellDepth; 1489a5b13a2SMatthew G Knepley cellDim -= depth - cellDepth; 1499a5b13a2SMatthew G Knepley ierr = PetscMalloc2(depth+1,PetscInt,&pStart,depth+1,PetscInt,&pEnd);CHKERRQ(ierr); 1509a5b13a2SMatthew G Knepley for (d = depth-1; d >= faceDepth; --d) { 1519a5b13a2SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, d, &pStart[d+1], &pEnd[d+1]);CHKERRQ(ierr); 1529f074e33SMatthew G Knepley } 1539a5b13a2SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, -1, NULL, &pStart[faceDepth]);CHKERRQ(ierr); 1549a5b13a2SMatthew G Knepley pEnd[faceDepth] = pStart[faceDepth]; 1559a5b13a2SMatthew G Knepley for (d = faceDepth-1; d >= 0; --d) { 1569a5b13a2SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 1579f074e33SMatthew G Knepley } 1589a5b13a2SMatthew G Knepley if (pEnd[cellDepth] > pStart[cellDepth]) {ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr);} 1599a5b13a2SMatthew 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); 1609a5b13a2SMatthew G Knepley ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 1619a5b13a2SMatthew G Knepley ierr = PetscHashIJKLSetMultivalued(faceTable, PETSC_FALSE);CHKERRQ(ierr); 1629a5b13a2SMatthew G Knepley for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) { 1639f074e33SMatthew G Knepley const PetscInt *cellFaces; 1649a5b13a2SMatthew G Knepley PetscInt numCellFaces, faceSize, cf, f; 1659f074e33SMatthew G Knepley 1669a5b13a2SMatthew G Knepley ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 1679a5b13a2SMatthew G Knepley if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 1689f074e33SMatthew G Knepley for (cf = 0; cf < numCellFaces; ++cf) { 1699a5b13a2SMatthew G Knepley const PetscInt *cellFace = &cellFaces[cf*faceSize]; 1709a5b13a2SMatthew G Knepley PetscHashIJKLKey key; 1719f074e33SMatthew G Knepley 1729a5b13a2SMatthew G Knepley if (faceSize == 2) { 1739a5b13a2SMatthew G Knepley key.i = PetscMin(cellFace[0], cellFace[1]); 1749a5b13a2SMatthew G Knepley key.j = PetscMax(cellFace[0], cellFace[1]); 1757fcd4ef0SJed Brown key.k = 0; 1767fcd4ef0SJed Brown key.l = 0; 1779a5b13a2SMatthew G Knepley } else { 1789a5b13a2SMatthew G Knepley key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0; 1799a5b13a2SMatthew G Knepley ierr = PetscSortInt(faceSize, (PetscInt *) &key); 1809f074e33SMatthew G Knepley } 1819a5b13a2SMatthew G Knepley ierr = PetscHashIJKLGet(faceTable, key, &f);CHKERRQ(ierr); 1829a5b13a2SMatthew G Knepley if (f < 0) { 1839a5b13a2SMatthew G Knepley ierr = PetscHashIJKLAdd(faceTable, key, face);CHKERRQ(ierr); 1849a5b13a2SMatthew G Knepley f = face++; 1859a5b13a2SMatthew G Knepley } 1869a5b13a2SMatthew G Knepley } 187439ece16SMatthew G. Knepley ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 1889a5b13a2SMatthew G Knepley } 1899a5b13a2SMatthew G Knepley pEnd[faceDepth] = face; 1909a5b13a2SMatthew G Knepley ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 1919a5b13a2SMatthew G Knepley /* Count new points */ 1929a5b13a2SMatthew G Knepley for (d = 0; d <= depth; ++d) { 1939a5b13a2SMatthew G Knepley numPoints += pEnd[d]-pStart[d]; 1949a5b13a2SMatthew G Knepley } 1959a5b13a2SMatthew G Knepley ierr = DMPlexSetChart(idm, 0, numPoints);CHKERRQ(ierr); 1969a5b13a2SMatthew G Knepley /* Set cone sizes */ 1979a5b13a2SMatthew G Knepley for (d = 0; d <= depth; ++d) { 1989a5b13a2SMatthew G Knepley PetscInt coneSize, p; 1999f074e33SMatthew G Knepley 2009a5b13a2SMatthew G Knepley if (d == faceDepth) { 2019a5b13a2SMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 2029a5b13a2SMatthew G Knepley /* I see no way to do this if we admit faces of different shapes */ 2039a5b13a2SMatthew G Knepley ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr); 2049a5b13a2SMatthew G Knepley } 205a014044eSMatthew G Knepley } else if (d == cellDepth) { 206a014044eSMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 207a014044eSMatthew G Knepley /* Number of cell faces may be different from number of cell vertices*/ 208a014044eSMatthew G Knepley ierr = DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);CHKERRQ(ierr); 209a014044eSMatthew G Knepley ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 210a014044eSMatthew G Knepley } 2119a5b13a2SMatthew G Knepley } else { 2129a5b13a2SMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 2139a5b13a2SMatthew G Knepley ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 2149a5b13a2SMatthew G Knepley ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 2159f074e33SMatthew G Knepley } 2169f074e33SMatthew G Knepley } 2179f074e33SMatthew G Knepley } 2189f074e33SMatthew G Knepley ierr = DMSetUp(idm);CHKERRQ(ierr); 2199f074e33SMatthew G Knepley /* Get face cones from subsets of cell vertices */ 2209a5b13a2SMatthew 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); 2219a5b13a2SMatthew G Knepley ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 2229a5b13a2SMatthew G Knepley ierr = PetscHashIJKLSetMultivalued(faceTable, PETSC_FALSE);CHKERRQ(ierr); 2239a5b13a2SMatthew G Knepley for (d = depth; d > cellDepth; --d) { 2249a5b13a2SMatthew G Knepley const PetscInt *cone; 2259a5b13a2SMatthew G Knepley PetscInt p; 2269a5b13a2SMatthew G Knepley 2279a5b13a2SMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 2289a5b13a2SMatthew G Knepley ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 2299a5b13a2SMatthew G Knepley ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr); 2309a5b13a2SMatthew G Knepley ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr); 2319a5b13a2SMatthew G Knepley ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr); 2329f074e33SMatthew G Knepley } 2339a5b13a2SMatthew G Knepley } 2349a5b13a2SMatthew G Knepley for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) { 2359f074e33SMatthew G Knepley const PetscInt *cellFaces; 2369a5b13a2SMatthew G Knepley PetscInt numCellFaces, faceSize, cf, f; 2379f074e33SMatthew G Knepley 2389a5b13a2SMatthew G Knepley ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 2399a5b13a2SMatthew G Knepley if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 2409f074e33SMatthew G Knepley for (cf = 0; cf < numCellFaces; ++cf) { 2419a5b13a2SMatthew G Knepley const PetscInt *cellFace = &cellFaces[cf*faceSize]; 2429a5b13a2SMatthew G Knepley PetscHashIJKLKey key; 2439f074e33SMatthew G Knepley 2449a5b13a2SMatthew G Knepley if (faceSize == 2) { 2459a5b13a2SMatthew G Knepley key.i = PetscMin(cellFace[0], cellFace[1]); 2469a5b13a2SMatthew G Knepley key.j = PetscMax(cellFace[0], cellFace[1]); 2477fcd4ef0SJed Brown key.k = 0; 2487fcd4ef0SJed Brown key.l = 0; 2499a5b13a2SMatthew G Knepley } else { 2509a5b13a2SMatthew G Knepley key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0; 2519a5b13a2SMatthew G Knepley ierr = PetscSortInt(faceSize, (PetscInt *) &key); 2529f074e33SMatthew G Knepley } 2539a5b13a2SMatthew G Knepley ierr = PetscHashIJKLGet(faceTable, key, &f);CHKERRQ(ierr); 2549a5b13a2SMatthew G Knepley if (f < 0) { 2559a5b13a2SMatthew G Knepley ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr); 2569a5b13a2SMatthew G Knepley ierr = PetscHashIJKLAdd(faceTable, key, face);CHKERRQ(ierr); 2579a5b13a2SMatthew G Knepley f = face++; 2589f074e33SMatthew G Knepley ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr); 2599a5b13a2SMatthew G Knepley } else { 2609a5b13a2SMatthew G Knepley const PetscInt *cone; 2619a5b13a2SMatthew G Knepley PetscInt coneSize, ornt, i, j; 2629f074e33SMatthew G Knepley 2639a5b13a2SMatthew G Knepley ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr); 2642e1b13c2SMatthew G. Knepley /* Orient face: Do not allow reverse orientation at the first vertex */ 2659f074e33SMatthew G Knepley ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr); 2669f074e33SMatthew G Knepley ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr); 2679a5b13a2SMatthew 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); 2689a5b13a2SMatthew G Knepley /* - First find the initial vertex */ 2699a5b13a2SMatthew G Knepley for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break; 2709a5b13a2SMatthew G Knepley /* - Try forward comparison */ 2719a5b13a2SMatthew G Knepley for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break; 2729a5b13a2SMatthew G Knepley if (j == faceSize) { 2739a5b13a2SMatthew G Knepley if ((faceSize == 2) && (i == 1)) ornt = -2; 2749a5b13a2SMatthew G Knepley else ornt = i; 2759a5b13a2SMatthew G Knepley } else { 2769a5b13a2SMatthew G Knepley /* - Try backward comparison */ 2779a5b13a2SMatthew G Knepley for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break; 2782e1b13c2SMatthew G. Knepley if (j == faceSize) { 2792e1b13c2SMatthew G. Knepley if (i == 0) ornt = -faceSize; 280dc1a705cSMatthew G. Knepley else ornt = -i; 2812e1b13c2SMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation"); 2829a5b13a2SMatthew G Knepley } 2839a5b13a2SMatthew G Knepley ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr); 2849f074e33SMatthew G Knepley } 2859f074e33SMatthew G Knepley } 286439ece16SMatthew G. Knepley ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 2879f074e33SMatthew G Knepley } 2889a5b13a2SMatthew 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]); 289c907b753SJed Brown ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 2909a5b13a2SMatthew G Knepley ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 2916551a8c7SMatthew G. Knepley ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 2929a5b13a2SMatthew G Knepley ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr); 2939a5b13a2SMatthew G Knepley ierr = DMPlexStratify(idm);CHKERRQ(ierr); 2949f074e33SMatthew G Knepley PetscFunctionReturn(0); 2959f074e33SMatthew G Knepley } 2969f074e33SMatthew G Knepley 2979f074e33SMatthew G Knepley #undef __FUNCT__ 2989f074e33SMatthew G Knepley #define __FUNCT__ "DMPlexInterpolate" 29980330477SMatthew G. Knepley /*@ 30080330477SMatthew G. Knepley DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc. 30180330477SMatthew G. Knepley 30280330477SMatthew G. Knepley Collective on DM 30380330477SMatthew G. Knepley 30480330477SMatthew G. Knepley Input Parameter: 3054e3744c5SMatthew G. Knepley . dm - The DMPlex object with only cells and vertices 30680330477SMatthew G. Knepley 30780330477SMatthew G. Knepley Output Parameter: 3084e3744c5SMatthew G. Knepley . dmInt - The complete DMPlex object 30980330477SMatthew G. Knepley 31080330477SMatthew G. Knepley Level: intermediate 31180330477SMatthew G. Knepley 31280330477SMatthew G. Knepley .keywords: mesh 3134e3744c5SMatthew G. Knepley .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList() 31480330477SMatthew G. Knepley @*/ 3159f074e33SMatthew G Knepley PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt) 3169f074e33SMatthew G Knepley { 3179a5b13a2SMatthew G Knepley DM idm, odm = dm; 3182e1b13c2SMatthew G. Knepley PetscInt depth, dim, d; 3199f074e33SMatthew G Knepley PetscErrorCode ierr; 3209f074e33SMatthew G Knepley 3219f074e33SMatthew G Knepley PetscFunctionBegin; 3222e1b13c2SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 3239f074e33SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 32476b791aaSMatthew G. Knepley if (dim <= 1) { 32576b791aaSMatthew G. Knepley ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 32676b791aaSMatthew G. Knepley idm = dm; 32776b791aaSMatthew G. Knepley } 3289a5b13a2SMatthew G Knepley for (d = 1; d < dim; ++d) { 3299a5b13a2SMatthew G Knepley /* Create interpolated mesh */ 3309a5b13a2SMatthew G Knepley ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr); 3319a5b13a2SMatthew G Knepley ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr); 3329a5b13a2SMatthew G Knepley ierr = DMPlexSetDimension(idm, dim);CHKERRQ(ierr); 3332e1b13c2SMatthew G. Knepley if (depth > 0) {ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr);} 3349a5b13a2SMatthew G Knepley if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);} 3359a5b13a2SMatthew G Knepley odm = idm; 3369f074e33SMatthew G Knepley } 3379a5b13a2SMatthew G Knepley *dmInt = idm; 3389f074e33SMatthew G Knepley PetscFunctionReturn(0); 3399f074e33SMatthew G Knepley } 34007b0a1fcSMatthew G Knepley 34107b0a1fcSMatthew G Knepley #undef __FUNCT__ 34207b0a1fcSMatthew G Knepley #define __FUNCT__ "DMPlexCopyCoordinates" 34380330477SMatthew G. Knepley /*@ 34480330477SMatthew G. Knepley DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices 34580330477SMatthew G. Knepley 34680330477SMatthew G. Knepley Collective on DM 34780330477SMatthew G. Knepley 34880330477SMatthew G. Knepley Input Parameter: 34980330477SMatthew G. Knepley . dmA - The DMPlex object with initial coordinates 35080330477SMatthew G. Knepley 35180330477SMatthew G. Knepley Output Parameter: 35280330477SMatthew G. Knepley . dmB - The DMPlex object with copied coordinates 35380330477SMatthew G. Knepley 35480330477SMatthew G. Knepley Level: intermediate 35580330477SMatthew G. Knepley 35680330477SMatthew G. Knepley Note: This is typically used when adding pieces other than vertices to a mesh 35780330477SMatthew G. Knepley 35880330477SMatthew G. Knepley .keywords: mesh 3595c386225SMatthew G. Knepley .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMPlexGetCoordinateSection() 36080330477SMatthew G. Knepley @*/ 36107b0a1fcSMatthew G Knepley PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB) 36207b0a1fcSMatthew G Knepley { 36307b0a1fcSMatthew G Knepley Vec coordinatesA, coordinatesB; 36407b0a1fcSMatthew G Knepley PetscSection coordSectionA, coordSectionB; 36507b0a1fcSMatthew G Knepley PetscScalar *coordsA, *coordsB; 36607b0a1fcSMatthew G Knepley PetscInt spaceDim, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d; 36707b0a1fcSMatthew G Knepley PetscErrorCode ierr; 36807b0a1fcSMatthew G Knepley 36907b0a1fcSMatthew G Knepley PetscFunctionBegin; 37076b791aaSMatthew G. Knepley if (dmA == dmB) PetscFunctionReturn(0); 37107b0a1fcSMatthew G Knepley ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr); 37207b0a1fcSMatthew G Knepley ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr); 37307b0a1fcSMatthew 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); 37407b0a1fcSMatthew G Knepley ierr = DMPlexGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr); 37507b0a1fcSMatthew G Knepley ierr = DMPlexGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr); 37607b0a1fcSMatthew G Knepley ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr); 37707b0a1fcSMatthew G Knepley ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr); 37807b0a1fcSMatthew G Knepley ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr); 37907b0a1fcSMatthew G Knepley ierr = PetscSectionSetChart(coordSectionB, vStartB, vEndB);CHKERRQ(ierr); 38007b0a1fcSMatthew G Knepley for (v = vStartB; v < vEndB; ++v) { 38107b0a1fcSMatthew G Knepley ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr); 38207b0a1fcSMatthew G Knepley ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr); 38307b0a1fcSMatthew G Knepley } 38407b0a1fcSMatthew G Knepley ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr); 38507b0a1fcSMatthew G Knepley ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr); 38607b0a1fcSMatthew G Knepley ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr); 38707b0a1fcSMatthew G Knepley ierr = VecCreate(PetscObjectComm((PetscObject) dmB), &coordinatesB);CHKERRQ(ierr); 38807b0a1fcSMatthew G Knepley ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr); 38907b0a1fcSMatthew G Knepley ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr); 3902eb5907fSJed Brown ierr = VecSetType(coordinatesB,VECSTANDARD);CHKERRQ(ierr); 39107b0a1fcSMatthew G Knepley ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr); 39207b0a1fcSMatthew G Knepley ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr); 39307b0a1fcSMatthew G Knepley for (v = 0; v < vEndB-vStartB; ++v) { 39407b0a1fcSMatthew G Knepley for (d = 0; d < spaceDim; ++d) { 39507b0a1fcSMatthew G Knepley coordsB[v*spaceDim+d] = coordsA[v*spaceDim+d]; 39607b0a1fcSMatthew G Knepley } 39707b0a1fcSMatthew G Knepley } 39807b0a1fcSMatthew G Knepley ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr); 39907b0a1fcSMatthew G Knepley ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr); 40007b0a1fcSMatthew G Knepley ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr); 40107b0a1fcSMatthew G Knepley ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr); 40207b0a1fcSMatthew G Knepley PetscFunctionReturn(0); 40307b0a1fcSMatthew G Knepley } 4045c386225SMatthew G. Knepley 4055c386225SMatthew G. Knepley #undef __FUNCT__ 4065c386225SMatthew G. Knepley #define __FUNCT__ "DMPlexCopyLabels" 4075c386225SMatthew G. Knepley /*@ 4085c386225SMatthew G. Knepley DMPlexCopyLabels - Copy labels from one mesh to another with a superset of the points 4095c386225SMatthew G. Knepley 4105c386225SMatthew G. Knepley Collective on DM 4115c386225SMatthew G. Knepley 4125c386225SMatthew G. Knepley Input Parameter: 4135c386225SMatthew G. Knepley . dmA - The DMPlex object with initial labels 4145c386225SMatthew G. Knepley 4155c386225SMatthew G. Knepley Output Parameter: 4165c386225SMatthew G. Knepley . dmB - The DMPlex object with copied labels 4175c386225SMatthew G. Knepley 4185c386225SMatthew G. Knepley Level: intermediate 4195c386225SMatthew G. Knepley 4205c386225SMatthew G. Knepley Note: This is typically used when interpolating or otherwise adding to a mesh 4215c386225SMatthew G. Knepley 4225c386225SMatthew G. Knepley .keywords: mesh 4235c386225SMatthew G. Knepley .seealso: DMCopyCoordinates(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMPlexGetCoordinateSection() 4245c386225SMatthew G. Knepley @*/ 4255c386225SMatthew G. Knepley PetscErrorCode DMPlexCopyLabels(DM dmA, DM dmB) 4265c386225SMatthew G. Knepley { 4275c386225SMatthew G. Knepley PetscInt numLabels, l; 4285c386225SMatthew G. Knepley PetscErrorCode ierr; 4295c386225SMatthew G. Knepley 4305c386225SMatthew G. Knepley PetscFunctionBegin; 4315c386225SMatthew G. Knepley if (dmA == dmB) PetscFunctionReturn(0); 4325c386225SMatthew G. Knepley ierr = DMPlexGetNumLabels(dmA, &numLabels);CHKERRQ(ierr); 4335c386225SMatthew G. Knepley for (l = 0; l < numLabels; ++l) { 4345c386225SMatthew G. Knepley DMLabel label, labelNew; 4355c386225SMatthew G. Knepley const char *name; 4365c386225SMatthew G. Knepley PetscBool flg; 4375c386225SMatthew G. Knepley 4385c386225SMatthew G. Knepley ierr = DMPlexGetLabelName(dmA, l, &name);CHKERRQ(ierr); 4395c386225SMatthew G. Knepley ierr = PetscStrcmp(name, "depth", &flg);CHKERRQ(ierr); 4405c386225SMatthew G. Knepley if (flg) continue; 4415c386225SMatthew G. Knepley ierr = DMPlexGetLabel(dmA, name, &label);CHKERRQ(ierr); 4425c386225SMatthew G. Knepley ierr = DMLabelDuplicate(label, &labelNew);CHKERRQ(ierr); 4435c386225SMatthew G. Knepley ierr = DMPlexAddLabel(dmB, labelNew);CHKERRQ(ierr); 4445c386225SMatthew G. Knepley } 4455c386225SMatthew G. Knepley PetscFunctionReturn(0); 4465c386225SMatthew G. Knepley } 447*f5cedc29SMatthew G. Knepley 448*f5cedc29SMatthew G. Knepley #undef __FUNCT__ 4494e3744c5SMatthew G. Knepley #define __FUNCT__ "DMPlexUninterpolate" 4504e3744c5SMatthew G. Knepley /*@ 4514e3744c5SMatthew G. Knepley DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh 4524e3744c5SMatthew G. Knepley 4534e3744c5SMatthew G. Knepley Collective on DM 4544e3744c5SMatthew G. Knepley 4554e3744c5SMatthew G. Knepley Input Parameter: 4564e3744c5SMatthew G. Knepley . dm - The complete DMPlex object 4574e3744c5SMatthew G. Knepley 4584e3744c5SMatthew G. Knepley Output Parameter: 4594e3744c5SMatthew G. Knepley . dmUnint - The DMPlex object with only cells and vertices 4604e3744c5SMatthew G. Knepley 4614e3744c5SMatthew G. Knepley Level: intermediate 4624e3744c5SMatthew G. Knepley 4634e3744c5SMatthew G. Knepley .keywords: mesh 4644e3744c5SMatthew G. Knepley .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList() 4654e3744c5SMatthew G. Knepley @*/ 4664e3744c5SMatthew G. Knepley PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint) 4674e3744c5SMatthew G. Knepley { 4684e3744c5SMatthew G. Knepley DM udm; 469595d4abbSMatthew G. Knepley PetscInt dim, vStart, vEnd, cStart, cEnd, c, maxConeSize = 0, *cone; 4704e3744c5SMatthew G. Knepley PetscErrorCode ierr; 4714e3744c5SMatthew G. Knepley 4724e3744c5SMatthew G. Knepley PetscFunctionBegin; 4734e3744c5SMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 4744e3744c5SMatthew G. Knepley if (dim <= 1) { 4754e3744c5SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 476595d4abbSMatthew G. Knepley *dmUnint = dm; 477595d4abbSMatthew G. Knepley PetscFunctionReturn(0); 4784e3744c5SMatthew G. Knepley } 479595d4abbSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 480595d4abbSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 4814e3744c5SMatthew G. Knepley ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr); 4824e3744c5SMatthew G. Knepley ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr); 4834e3744c5SMatthew G. Knepley ierr = DMPlexSetDimension(udm, dim);CHKERRQ(ierr); 4844e3744c5SMatthew G. Knepley ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr); 4854e3744c5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 486595d4abbSMatthew G. Knepley PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 4874e3744c5SMatthew G. Knepley 4884e3744c5SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 4894e3744c5SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 4904e3744c5SMatthew G. Knepley const PetscInt p = closure[cl]; 4914e3744c5SMatthew G. Knepley 4924e3744c5SMatthew G. Knepley if ((p >= vStart) && (p < vEnd)) ++coneSize; 4934e3744c5SMatthew G. Knepley } 4944e3744c5SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 4954e3744c5SMatthew G. Knepley ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr); 496595d4abbSMatthew G. Knepley maxConeSize = PetscMax(maxConeSize, coneSize); 4974e3744c5SMatthew G. Knepley } 4984e3744c5SMatthew G. Knepley ierr = DMSetUp(udm);CHKERRQ(ierr); 4994e3744c5SMatthew G. Knepley ierr = PetscMalloc(maxConeSize * sizeof(PetscInt), &cone);CHKERRQ(ierr); 5004e3744c5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 501595d4abbSMatthew G. Knepley PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 5024e3744c5SMatthew G. Knepley 5034e3744c5SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 5044e3744c5SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 5054e3744c5SMatthew G. Knepley const PetscInt p = closure[cl]; 5064e3744c5SMatthew G. Knepley 5074e3744c5SMatthew G. Knepley if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p; 5084e3744c5SMatthew G. Knepley } 5094e3744c5SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 5104e3744c5SMatthew G. Knepley ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr); 5114e3744c5SMatthew G. Knepley } 5124e3744c5SMatthew G. Knepley ierr = PetscFree(cone);CHKERRQ(ierr); 5134e3744c5SMatthew G. Knepley ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr); 5144e3744c5SMatthew G. Knepley ierr = DMPlexStratify(udm);CHKERRQ(ierr); 5154e3744c5SMatthew G. Knepley *dmUnint = udm; 5164e3744c5SMatthew G. Knepley PetscFunctionReturn(0); 5174e3744c5SMatthew G. Knepley } 518