1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 29f074e33SMatthew G Knepley #include <../src/sys/utils/hash.h> 39f074e33SMatthew G Knepley 49f074e33SMatthew G Knepley /* 59a5b13a2SMatthew G Knepley DMPlexGetFaces_Internal - Gets groups of vertices that correspond to faces for the given cell 69f074e33SMatthew G Knepley */ 7439ece16SMatthew G. Knepley PetscErrorCode DMPlexGetFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 89f074e33SMatthew G Knepley { 99f074e33SMatthew G Knepley const PetscInt *cone = NULL; 109a5b13a2SMatthew G Knepley PetscInt maxConeSize, maxSupportSize, coneSize; 119f074e33SMatthew G Knepley PetscErrorCode ierr; 129f074e33SMatthew G Knepley 139f074e33SMatthew G Knepley PetscFunctionBegin; 149f074e33SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 159a5b13a2SMatthew G Knepley ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 169f074e33SMatthew G Knepley ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 179f074e33SMatthew G Knepley ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 18439ece16SMatthew G. Knepley ierr = DMPlexGetRawFaces_Internal(dm, dim, coneSize, cone, numFaces, faceSize, faces);CHKERRQ(ierr); 19439ece16SMatthew G. Knepley PetscFunctionReturn(0); 20439ece16SMatthew G. Knepley } 21439ece16SMatthew G. Knepley 22439ece16SMatthew G. Knepley /* 23439ece16SMatthew G. Knepley DMPlexRestoreFaces_Internal - Restores the array 24439ece16SMatthew G. Knepley */ 25439ece16SMatthew G. Knepley PetscErrorCode DMPlexRestoreFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 26439ece16SMatthew G. Knepley { 27439ece16SMatthew G. Knepley PetscErrorCode ierr; 28439ece16SMatthew G. Knepley 29439ece16SMatthew G. Knepley PetscFunctionBegin; 3006acb219SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, 0, PETSC_INT, (void *) faces);CHKERRQ(ierr); 31439ece16SMatthew G. Knepley PetscFunctionReturn(0); 32439ece16SMatthew G. Knepley } 33439ece16SMatthew G. Knepley 34439ece16SMatthew G. Knepley /* 35439ece16SMatthew G. Knepley DMPlexGetRawFaces_Internal - Gets groups of vertices that correspond to faces for the given cone 36439ece16SMatthew G. Knepley */ 37439ece16SMatthew G. Knepley PetscErrorCode DMPlexGetRawFaces_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 38439ece16SMatthew G. Knepley { 39439ece16SMatthew G. Knepley PetscInt *facesTmp; 40439ece16SMatthew G. Knepley PetscInt maxConeSize, maxSupportSize; 41439ece16SMatthew G. Knepley PetscErrorCode ierr; 42439ece16SMatthew G. Knepley 43439ece16SMatthew G. Knepley PetscFunctionBegin; 44439ece16SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 45439ece16SMatthew G. Knepley ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 468496a43fSMatthew G. Knepley if (faces) {ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), PETSC_INT, &facesTmp);CHKERRQ(ierr);} 479f074e33SMatthew G Knepley switch (dim) { 48c49d9212SMatthew G. Knepley case 1: 49c49d9212SMatthew G. Knepley switch (coneSize) { 50c49d9212SMatthew G. Knepley case 2: 51c49d9212SMatthew G. Knepley if (faces) { 52c49d9212SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 53c49d9212SMatthew G. Knepley *faces = facesTmp; 54c49d9212SMatthew G. Knepley } 55c49d9212SMatthew G. Knepley if (numFaces) *numFaces = 2; 56c49d9212SMatthew G. Knepley if (faceSize) *faceSize = 1; 57c49d9212SMatthew G. Knepley break; 58c49d9212SMatthew G. Knepley default: 59c49d9212SMatthew G. Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 60c49d9212SMatthew G. Knepley } 61c49d9212SMatthew G. Knepley break; 629f074e33SMatthew G Knepley case 2: 639f074e33SMatthew G Knepley switch (coneSize) { 649f074e33SMatthew G Knepley case 3: 659a5b13a2SMatthew G Knepley if (faces) { 669a5b13a2SMatthew G Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 679a5b13a2SMatthew G Knepley facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 689a5b13a2SMatthew G Knepley facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 699a5b13a2SMatthew G Knepley *faces = facesTmp; 709a5b13a2SMatthew G Knepley } 719a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 3; 729a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 2; 739f074e33SMatthew G Knepley break; 749f074e33SMatthew G Knepley case 4: 759a5b13a2SMatthew G Knepley /* Vertices follow right hand rule */ 769a5b13a2SMatthew G Knepley if (faces) { 779a5b13a2SMatthew G Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 789a5b13a2SMatthew G Knepley facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 799a5b13a2SMatthew G Knepley facesTmp[4] = cone[2]; facesTmp[5] = cone[3]; 809a5b13a2SMatthew G Knepley facesTmp[6] = cone[3]; facesTmp[7] = cone[0]; 819a5b13a2SMatthew G Knepley *faces = facesTmp; 829a5b13a2SMatthew G Knepley } 839a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 4; 849a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 2; 859f074e33SMatthew G Knepley break; 869f074e33SMatthew G Knepley default: 879f074e33SMatthew G Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 889f074e33SMatthew G Knepley } 899f074e33SMatthew G Knepley break; 909f074e33SMatthew G Knepley case 3: 919f074e33SMatthew G Knepley switch (coneSize) { 929f074e33SMatthew G Knepley case 3: 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[0]; 979a5b13a2SMatthew G Knepley *faces = facesTmp; 989a5b13a2SMatthew G Knepley } 999a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 3; 1009a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 2; 1019f074e33SMatthew G Knepley break; 1029f074e33SMatthew G Knepley case 4: 1032e1b13c2SMatthew G. Knepley /* Vertices of first face follow right hand rule and normal points away from last vertex */ 1049a5b13a2SMatthew G Knepley if (faces) { 1052e1b13c2SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; 1062e1b13c2SMatthew G. Knepley facesTmp[3] = cone[0]; facesTmp[4] = cone[3]; facesTmp[5] = cone[1]; 1072e1b13c2SMatthew G. Knepley facesTmp[6] = cone[0]; facesTmp[7] = cone[2]; facesTmp[8] = cone[3]; 1082e1b13c2SMatthew G. Knepley facesTmp[9] = cone[2]; facesTmp[10] = cone[1]; facesTmp[11] = cone[3]; 1099a5b13a2SMatthew G Knepley *faces = facesTmp; 1109a5b13a2SMatthew G Knepley } 1119a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 4; 1129a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 3; 1139a5b13a2SMatthew G Knepley break; 1149a5b13a2SMatthew G Knepley case 8: 1159a5b13a2SMatthew G Knepley if (faces) { 11651cfd6a4SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; facesTmp[3] = cone[3]; /* Bottom */ 11751cfd6a4SMatthew G. Knepley facesTmp[4] = cone[4]; facesTmp[5] = cone[5]; facesTmp[6] = cone[6]; facesTmp[7] = cone[7]; /* Top */ 11851cfd6a4SMatthew G. Knepley facesTmp[8] = cone[0]; facesTmp[9] = cone[3]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4]; /* Front */ 11951cfd6a4SMatthew G. Knepley facesTmp[12] = cone[2]; facesTmp[13] = cone[1]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6]; /* Back */ 12051cfd6a4SMatthew G. Knepley facesTmp[16] = cone[3]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5]; /* Right */ 12151cfd6a4SMatthew G. Knepley facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[1]; /* Left */ 1229a5b13a2SMatthew G Knepley *faces = facesTmp; 1239a5b13a2SMatthew G Knepley } 1249a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 6; 1259a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 4; 1269f074e33SMatthew G Knepley break; 1279f074e33SMatthew G Knepley default: 1289f074e33SMatthew G Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 1299f074e33SMatthew G Knepley } 1309f074e33SMatthew G Knepley break; 1319f074e33SMatthew G Knepley default: 1329f074e33SMatthew G Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim); 1339f074e33SMatthew G Knepley } 1349f074e33SMatthew G Knepley PetscFunctionReturn(0); 1359f074e33SMatthew G Knepley } 1369f074e33SMatthew G Knepley 1379a5b13a2SMatthew G Knepley /* This interpolates faces for cells at some stratum */ 1389a5b13a2SMatthew G Knepley static PetscErrorCode DMPlexInterpolateFaces_Internal(DM dm, PetscInt cellDepth, DM idm) 1399f074e33SMatthew G Knepley { 140d4efc6f1SMatthew G. Knepley DMLabel subpointMap; 1419a5b13a2SMatthew G Knepley PetscHashIJKL faceTable; 1429a5b13a2SMatthew G Knepley PetscInt *pStart, *pEnd; 1439a5b13a2SMatthew G Knepley PetscInt cellDim, depth, faceDepth = cellDepth, numPoints = 0, faceSizeAll = 0, face, c, d; 1449f074e33SMatthew G Knepley PetscErrorCode ierr; 1459f074e33SMatthew G Knepley 1469f074e33SMatthew G Knepley PetscFunctionBegin; 147c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &cellDim);CHKERRQ(ierr); 148d4efc6f1SMatthew G. Knepley /* HACK: I need a better way to determine face dimension, or an alternative to GetFaces() */ 149d4efc6f1SMatthew G. Knepley ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr); 150d4efc6f1SMatthew G. Knepley if (subpointMap) ++cellDim; 1519a5b13a2SMatthew G Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1529a5b13a2SMatthew G Knepley ++depth; 1539a5b13a2SMatthew G Knepley ++cellDepth; 1549a5b13a2SMatthew G Knepley cellDim -= depth - cellDepth; 155dcca6d9dSJed Brown ierr = PetscMalloc2(depth+1,&pStart,depth+1,&pEnd);CHKERRQ(ierr); 1569a5b13a2SMatthew G Knepley for (d = depth-1; d >= faceDepth; --d) { 1579a5b13a2SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, d, &pStart[d+1], &pEnd[d+1]);CHKERRQ(ierr); 1589f074e33SMatthew G Knepley } 1599a5b13a2SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, -1, NULL, &pStart[faceDepth]);CHKERRQ(ierr); 1609a5b13a2SMatthew G Knepley pEnd[faceDepth] = pStart[faceDepth]; 1619a5b13a2SMatthew G Knepley for (d = faceDepth-1; d >= 0; --d) { 1629a5b13a2SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 1639f074e33SMatthew G Knepley } 1649a5b13a2SMatthew G Knepley if (pEnd[cellDepth] > pStart[cellDepth]) {ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr);} 1659a5b13a2SMatthew 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); 1669a5b13a2SMatthew G Knepley ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 1679a5b13a2SMatthew G Knepley for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) { 1689f074e33SMatthew G Knepley const PetscInt *cellFaces; 169735a0255SMatthew G. Knepley PetscInt numCellFaces, faceSize, cf; 1709f074e33SMatthew G Knepley 1719a5b13a2SMatthew G Knepley ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 1729a5b13a2SMatthew G Knepley if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 1739f074e33SMatthew G Knepley for (cf = 0; cf < numCellFaces; ++cf) { 1749a5b13a2SMatthew G Knepley const PetscInt *cellFace = &cellFaces[cf*faceSize]; 1759a5b13a2SMatthew G Knepley PetscHashIJKLKey key; 176735a0255SMatthew G. Knepley PetscHashIJKLIter missing, iter; 1779f074e33SMatthew G Knepley 1789a5b13a2SMatthew G Knepley if (faceSize == 2) { 1799a5b13a2SMatthew G Knepley key.i = PetscMin(cellFace[0], cellFace[1]); 1809a5b13a2SMatthew G Knepley key.j = PetscMax(cellFace[0], cellFace[1]); 1817fcd4ef0SJed Brown key.k = 0; 1827fcd4ef0SJed Brown key.l = 0; 1839a5b13a2SMatthew G Knepley } else { 1849a5b13a2SMatthew G Knepley key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0; 185302440fdSBarry Smith ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr); 1869f074e33SMatthew G Knepley } 187735a0255SMatthew G. Knepley ierr = PetscHashIJKLPut(faceTable, key, &missing, &iter);CHKERRQ(ierr); 188735a0255SMatthew G. Knepley if (missing) {ierr = PetscHashIJKLSet(faceTable, iter, face++);CHKERRQ(ierr);} 1899a5b13a2SMatthew G Knepley } 190439ece16SMatthew G. Knepley ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 1919a5b13a2SMatthew G Knepley } 1929a5b13a2SMatthew G Knepley pEnd[faceDepth] = face; 1939a5b13a2SMatthew G Knepley ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 1949a5b13a2SMatthew G Knepley /* Count new points */ 1959a5b13a2SMatthew G Knepley for (d = 0; d <= depth; ++d) { 1969a5b13a2SMatthew G Knepley numPoints += pEnd[d]-pStart[d]; 1979a5b13a2SMatthew G Knepley } 1989a5b13a2SMatthew G Knepley ierr = DMPlexSetChart(idm, 0, numPoints);CHKERRQ(ierr); 1999a5b13a2SMatthew G Knepley /* Set cone sizes */ 2009a5b13a2SMatthew G Knepley for (d = 0; d <= depth; ++d) { 2019a5b13a2SMatthew G Knepley PetscInt coneSize, p; 2029f074e33SMatthew G Knepley 2039a5b13a2SMatthew G Knepley if (d == faceDepth) { 2049a5b13a2SMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 2059a5b13a2SMatthew G Knepley /* I see no way to do this if we admit faces of different shapes */ 2069a5b13a2SMatthew G Knepley ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr); 2079a5b13a2SMatthew G Knepley } 208a014044eSMatthew G Knepley } else if (d == cellDepth) { 209a014044eSMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 210a014044eSMatthew G Knepley /* Number of cell faces may be different from number of cell vertices*/ 211a014044eSMatthew G Knepley ierr = DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);CHKERRQ(ierr); 212a014044eSMatthew G Knepley ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 213a014044eSMatthew G Knepley } 2149a5b13a2SMatthew G Knepley } else { 2159a5b13a2SMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 2169a5b13a2SMatthew G Knepley ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 2179a5b13a2SMatthew G Knepley ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 2189f074e33SMatthew G Knepley } 2199f074e33SMatthew G Knepley } 2209f074e33SMatthew G Knepley } 2219f074e33SMatthew G Knepley ierr = DMSetUp(idm);CHKERRQ(ierr); 2229f074e33SMatthew G Knepley /* Get face cones from subsets of cell vertices */ 2239a5b13a2SMatthew 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); 2249a5b13a2SMatthew G Knepley ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 2259a5b13a2SMatthew G Knepley for (d = depth; d > cellDepth; --d) { 2269a5b13a2SMatthew G Knepley const PetscInt *cone; 2279a5b13a2SMatthew G Knepley PetscInt p; 2289a5b13a2SMatthew G Knepley 2299a5b13a2SMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 2309a5b13a2SMatthew G Knepley ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 2319a5b13a2SMatthew G Knepley ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr); 2329a5b13a2SMatthew G Knepley ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr); 2339a5b13a2SMatthew G Knepley ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr); 2349f074e33SMatthew G Knepley } 2359a5b13a2SMatthew G Knepley } 2369a5b13a2SMatthew G Knepley for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) { 2379f074e33SMatthew G Knepley const PetscInt *cellFaces; 238735a0255SMatthew G. Knepley PetscInt numCellFaces, faceSize, cf; 2399f074e33SMatthew G Knepley 2409a5b13a2SMatthew G Knepley ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 2419a5b13a2SMatthew G Knepley if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 2429f074e33SMatthew G Knepley for (cf = 0; cf < numCellFaces; ++cf) { 2439a5b13a2SMatthew G Knepley const PetscInt *cellFace = &cellFaces[cf*faceSize]; 2449a5b13a2SMatthew G Knepley PetscHashIJKLKey key; 245735a0255SMatthew G. Knepley PetscHashIJKLIter missing, iter; 2469f074e33SMatthew G Knepley 2479a5b13a2SMatthew G Knepley if (faceSize == 2) { 2489a5b13a2SMatthew G Knepley key.i = PetscMin(cellFace[0], cellFace[1]); 2499a5b13a2SMatthew G Knepley key.j = PetscMax(cellFace[0], cellFace[1]); 2507fcd4ef0SJed Brown key.k = 0; 2517fcd4ef0SJed Brown key.l = 0; 2529a5b13a2SMatthew G Knepley } else { 2539a5b13a2SMatthew G Knepley key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0; 254302440fdSBarry Smith ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr); 2559f074e33SMatthew G Knepley } 256735a0255SMatthew G. Knepley ierr = PetscHashIJKLPut(faceTable, key, &missing, &iter);CHKERRQ(ierr); 257735a0255SMatthew G. Knepley if (missing) { 2589a5b13a2SMatthew G Knepley ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr); 259735a0255SMatthew G. Knepley ierr = PetscHashIJKLSet(faceTable, iter, face);CHKERRQ(ierr); 260735a0255SMatthew G. Knepley ierr = DMPlexInsertCone(idm, c, cf, face++);CHKERRQ(ierr); 2619a5b13a2SMatthew G Knepley } else { 2629a5b13a2SMatthew G Knepley const PetscInt *cone; 263735a0255SMatthew G. Knepley PetscInt coneSize, ornt, i, j, f; 2649f074e33SMatthew G Knepley 265735a0255SMatthew G. Knepley ierr = PetscHashIJKLGet(faceTable, iter, &f);CHKERRQ(ierr); 2669a5b13a2SMatthew G Knepley ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr); 2672e1b13c2SMatthew G. Knepley /* Orient face: Do not allow reverse orientation at the first vertex */ 2689f074e33SMatthew G Knepley ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr); 2699f074e33SMatthew G Knepley ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr); 2709a5b13a2SMatthew 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); 2719a5b13a2SMatthew G Knepley /* - First find the initial vertex */ 2729a5b13a2SMatthew G Knepley for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break; 2739a5b13a2SMatthew G Knepley /* - Try forward comparison */ 2749a5b13a2SMatthew G Knepley for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break; 2759a5b13a2SMatthew G Knepley if (j == faceSize) { 2769a5b13a2SMatthew G Knepley if ((faceSize == 2) && (i == 1)) ornt = -2; 2779a5b13a2SMatthew G Knepley else ornt = i; 2789a5b13a2SMatthew G Knepley } else { 2799a5b13a2SMatthew G Knepley /* - Try backward comparison */ 2809a5b13a2SMatthew G Knepley for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break; 2812e1b13c2SMatthew G. Knepley if (j == faceSize) { 2822e1b13c2SMatthew G. Knepley if (i == 0) ornt = -faceSize; 283dc1a705cSMatthew G. Knepley else ornt = -i; 2842e1b13c2SMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation"); 2859a5b13a2SMatthew G Knepley } 2869a5b13a2SMatthew G Knepley ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr); 2879f074e33SMatthew G Knepley } 2889f074e33SMatthew G Knepley } 289439ece16SMatthew G. Knepley ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 2909f074e33SMatthew G Knepley } 2919a5b13a2SMatthew 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]); 292c907b753SJed Brown ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 2939a5b13a2SMatthew G Knepley ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 2946551a8c7SMatthew G. Knepley ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 2959a5b13a2SMatthew G Knepley ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr); 2969a5b13a2SMatthew G Knepley ierr = DMPlexStratify(idm);CHKERRQ(ierr); 2979f074e33SMatthew G Knepley PetscFunctionReturn(0); 2989f074e33SMatthew G Knepley } 2999f074e33SMatthew G Knepley 3007bffde78SMichael Lange /* This interpolates the PointSF in parallel following local interpolation */ 3017bffde78SMichael Lange static PetscErrorCode DMPlexInterpolatePointSF(DM dm, PetscSF pointSF, PetscInt depth) 3027bffde78SMichael Lange { 3039852e123SBarry Smith PetscMPIInt size, rank; 3047bffde78SMichael Lange PetscInt p, c, d, dof, offset; 305ffd6914dSSatish Balay PetscInt numLeaves, numRoots, candidatesSize, candidatesRemoteSize; 3067bffde78SMichael Lange const PetscInt *localPoints; 3077bffde78SMichael Lange const PetscSFNode *remotePoints; 3087bffde78SMichael Lange PetscSFNode *candidates, *candidatesRemote, *claims; 3097bffde78SMichael Lange PetscSection candidateSection, candidateSectionRemote, claimSection; 3107bffde78SMichael Lange PetscHashI leafhash; 3117bffde78SMichael Lange PetscHashIJ roothash; 3127bffde78SMichael Lange PetscHashIJKey key; 3137bffde78SMichael Lange PetscErrorCode ierr; 3147bffde78SMichael Lange 3157bffde78SMichael Lange PetscFunctionBegin; 3169852e123SBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr); 3177bffde78SMichael Lange ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 3187bffde78SMichael Lange ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 3199852e123SBarry Smith if (size < 2 || numRoots < 0) PetscFunctionReturn(0); 32025afeb17SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr); 3217bffde78SMichael Lange /* Build hashes of points in the SF for efficient lookup */ 3227bffde78SMichael Lange PetscHashICreate(leafhash); 3237bffde78SMichael Lange PetscHashIJCreate(&roothash); 3247bffde78SMichael Lange ierr = PetscHashIJSetMultivalued(roothash, PETSC_FALSE);CHKERRQ(ierr); 3257bffde78SMichael Lange for (p = 0; p < numLeaves; ++p) { 3267bffde78SMichael Lange PetscHashIAdd(leafhash, localPoints[p], p); 3277bffde78SMichael Lange key.i = remotePoints[p].index; key.j = remotePoints[p].rank; 3287bffde78SMichael Lange PetscHashIJAdd(roothash, key, p); 3297bffde78SMichael Lange } 3307bffde78SMichael Lange /* Build a section / SFNode array of candidate points in the single-level adjacency of leaves, 3317bffde78SMichael Lange where each candidate is defined by the root entry for the other vertex that defines the edge. */ 3327bffde78SMichael Lange ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSection);CHKERRQ(ierr); 3337bffde78SMichael Lange ierr = PetscSectionSetChart(candidateSection, 0, numRoots);CHKERRQ(ierr); 3347bffde78SMichael Lange { 335ffd6914dSSatish Balay PetscInt leaf, root, idx, a, *adj = NULL; 3367bffde78SMichael Lange for (p = 0; p < numLeaves; ++p) { 3377bffde78SMichael Lange PetscInt adjSize = PETSC_DETERMINE; 3387bffde78SMichael Lange ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr); 3397bffde78SMichael Lange for (a = 0; a < adjSize; ++a) { 3407bffde78SMichael Lange PetscHashIMap(leafhash, adj[a], leaf); 3417bffde78SMichael Lange if (leaf >= 0) {ierr = PetscSectionAddDof(candidateSection, localPoints[p], 1);CHKERRQ(ierr);} 3427bffde78SMichael Lange } 3437bffde78SMichael Lange } 3447bffde78SMichael Lange ierr = PetscSectionSetUp(candidateSection);CHKERRQ(ierr); 3457bffde78SMichael Lange ierr = PetscSectionGetStorageSize(candidateSection, &candidatesSize);CHKERRQ(ierr); 3467bffde78SMichael Lange ierr = PetscMalloc1(candidatesSize, &candidates);CHKERRQ(ierr); 3477bffde78SMichael Lange for (p = 0; p < numLeaves; ++p) { 3487bffde78SMichael Lange PetscInt adjSize = PETSC_DETERMINE; 3497bffde78SMichael Lange ierr = PetscSectionGetOffset(candidateSection, localPoints[p], &offset);CHKERRQ(ierr); 3507bffde78SMichael Lange ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr); 3517bffde78SMichael Lange for (idx = 0, a = 0; a < adjSize; ++a) { 3527bffde78SMichael Lange PetscHashIMap(leafhash, adj[a], root); 3537bffde78SMichael Lange if (root >= 0) candidates[offset+idx++] = remotePoints[root]; 3547bffde78SMichael Lange } 3557bffde78SMichael Lange } 3567bffde78SMichael Lange ierr = PetscFree(adj);CHKERRQ(ierr); 3577bffde78SMichael Lange } 3587bffde78SMichael Lange /* Gather candidate section / array pair into the root partition via inverse(multi(pointSF)). */ 3597bffde78SMichael Lange { 3607bffde78SMichael Lange PetscSF sfMulti, sfInverse, sfCandidates; 3617bffde78SMichael Lange PetscInt *remoteOffsets; 3627bffde78SMichael Lange ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr); 3637bffde78SMichael Lange ierr = PetscSFCreateInverseSF(sfMulti, &sfInverse);CHKERRQ(ierr); 3647bffde78SMichael Lange ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSectionRemote);CHKERRQ(ierr); 3657bffde78SMichael Lange ierr = PetscSFDistributeSection(sfInverse, candidateSection, &remoteOffsets, candidateSectionRemote);CHKERRQ(ierr); 3667bffde78SMichael Lange ierr = PetscSFCreateSectionSF(sfInverse, candidateSection, remoteOffsets, candidateSectionRemote, &sfCandidates);CHKERRQ(ierr); 3677bffde78SMichael Lange ierr = PetscSectionGetStorageSize(candidateSectionRemote, &candidatesRemoteSize);CHKERRQ(ierr); 3687bffde78SMichael Lange ierr = PetscMalloc1(candidatesRemoteSize, &candidatesRemote);CHKERRQ(ierr); 3697bffde78SMichael Lange ierr = PetscSFBcastBegin(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr); 3707bffde78SMichael Lange ierr = PetscSFBcastEnd(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr); 3717bffde78SMichael Lange ierr = PetscSFDestroy(&sfInverse);CHKERRQ(ierr); 3727bffde78SMichael Lange ierr = PetscSFDestroy(&sfCandidates);CHKERRQ(ierr); 3737bffde78SMichael Lange ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 3747bffde78SMichael Lange } 3757bffde78SMichael Lange /* Walk local roots and check for each remote candidate whether we know all required points, 3767bffde78SMichael Lange either from owning it or having a root entry in the point SF. If we do we place a claim 3777bffde78SMichael Lange by replacing the vertex number with our edge ID. */ 3787bffde78SMichael Lange { 3797bffde78SMichael Lange PetscInt idx, root, joinSize, vertices[2]; 3807bffde78SMichael Lange const PetscInt *rootdegree, *join = NULL; 3817bffde78SMichael Lange ierr = PetscSFComputeDegreeBegin(pointSF, &rootdegree);CHKERRQ(ierr); 3827bffde78SMichael Lange ierr = PetscSFComputeDegreeEnd(pointSF, &rootdegree);CHKERRQ(ierr); 3837bffde78SMichael Lange /* Loop remote edge connections and put in a claim if both vertices are known */ 3847bffde78SMichael Lange for (idx = 0, p = 0; p < numRoots; ++p) { 3857bffde78SMichael Lange for (d = 0; d < rootdegree[p]; ++d) { 3867bffde78SMichael Lange ierr = PetscSectionGetDof(candidateSectionRemote, idx, &dof);CHKERRQ(ierr); 3877bffde78SMichael Lange ierr = PetscSectionGetOffset(candidateSectionRemote, idx, &offset);CHKERRQ(ierr); 3887bffde78SMichael Lange for (c = 0; c < dof; ++c) { 3897bffde78SMichael Lange /* We own both vertices, so we claim the edge by replacing vertex with edge */ 3907bffde78SMichael Lange if (candidatesRemote[offset+c].rank == rank) { 3917bffde78SMichael Lange vertices[0] = p; vertices[1] = candidatesRemote[offset+c].index; 3927bffde78SMichael Lange ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 3937bffde78SMichael Lange if (joinSize == 1) candidatesRemote[offset+c].index = join[0]; 3947bffde78SMichael Lange ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 3957bffde78SMichael Lange continue; 3967bffde78SMichael Lange } 3977bffde78SMichael Lange /* If we own one vertex and share a root with the other, we claim it */ 3987bffde78SMichael Lange key.i = candidatesRemote[offset+c].index; key.j = candidatesRemote[offset+c].rank; 3997bffde78SMichael Lange PetscHashIJGet(roothash, key, &root); 4007bffde78SMichael Lange if (root >= 0) { 4017bffde78SMichael Lange vertices[0] = p; vertices[1] = localPoints[root]; 4027bffde78SMichael Lange ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 4037bffde78SMichael Lange if (joinSize == 1) { 4047bffde78SMichael Lange candidatesRemote[offset+c].index = join[0]; 4057bffde78SMichael Lange candidatesRemote[offset+c].rank = rank; 4067bffde78SMichael Lange } 4077bffde78SMichael Lange ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 4087bffde78SMichael Lange } 4097bffde78SMichael Lange } 4107bffde78SMichael Lange idx++; 4117bffde78SMichael Lange } 4127bffde78SMichael Lange } 4137bffde78SMichael Lange } 4147bffde78SMichael Lange /* Push claims back to receiver via the MultiSF and derive new pointSF mapping on receiver */ 4157bffde78SMichael Lange { 4167bffde78SMichael Lange PetscSF sfMulti, sfClaims, sfPointNew; 4177bffde78SMichael Lange PetscHashI claimshash; 4187bffde78SMichael Lange PetscInt size, pStart, pEnd, root, joinSize, numLocalNew; 4197bffde78SMichael Lange PetscInt *remoteOffsets, *localPointsNew, vertices[2]; 420ffd6914dSSatish Balay const PetscInt *join = NULL; 4217bffde78SMichael Lange PetscSFNode *remotePointsNew; 4227bffde78SMichael Lange ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr); 4237bffde78SMichael Lange ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &claimSection);CHKERRQ(ierr); 4247bffde78SMichael Lange ierr = PetscSFDistributeSection(sfMulti, candidateSectionRemote, &remoteOffsets, claimSection);CHKERRQ(ierr); 4257bffde78SMichael Lange ierr = PetscSFCreateSectionSF(sfMulti, candidateSectionRemote, remoteOffsets, claimSection, &sfClaims);CHKERRQ(ierr); 4267bffde78SMichael Lange ierr = PetscSectionGetStorageSize(claimSection, &size);CHKERRQ(ierr); 4277bffde78SMichael Lange ierr = PetscMalloc1(size, &claims);CHKERRQ(ierr); 4287bffde78SMichael Lange ierr = PetscSFBcastBegin(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr); 4297bffde78SMichael Lange ierr = PetscSFBcastEnd(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr); 4307bffde78SMichael Lange ierr = PetscSFDestroy(&sfClaims);CHKERRQ(ierr); 4317bffde78SMichael Lange ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 4327bffde78SMichael Lange /* Walk the original section of local supports and add an SF entry for each updated item */ 4337bffde78SMichael Lange PetscHashICreate(claimshash); 4347bffde78SMichael Lange for (p = 0; p < numRoots; ++p) { 4357bffde78SMichael Lange ierr = PetscSectionGetDof(candidateSection, p, &dof);CHKERRQ(ierr); 4367bffde78SMichael Lange ierr = PetscSectionGetOffset(candidateSection, p, &offset);CHKERRQ(ierr); 4377bffde78SMichael Lange for (d = 0; d < dof; ++d) { 4387bffde78SMichael Lange if (candidates[offset+d].index != claims[offset+d].index) { 4397bffde78SMichael Lange key.i = candidates[offset+d].index; key.j = candidates[offset+d].rank; 4407bffde78SMichael Lange PetscHashIJGet(roothash, key, &root); 4417bffde78SMichael Lange if (root >= 0) { 4427bffde78SMichael Lange vertices[0] = p; vertices[1] = localPoints[root]; 4437bffde78SMichael Lange ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 4447bffde78SMichael Lange if (joinSize == 1) PetscHashIAdd(claimshash, join[0], offset+d); 4457bffde78SMichael Lange ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 4467bffde78SMichael Lange } 4477bffde78SMichael Lange } 4487bffde78SMichael Lange } 4497bffde78SMichael Lange } 4507bffde78SMichael Lange /* Create new pointSF from hashed claims */ 4517bffde78SMichael Lange PetscHashISize(claimshash, numLocalNew); 4527bffde78SMichael Lange ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 4537bffde78SMichael Lange ierr = PetscMalloc1(numLeaves + numLocalNew, &localPointsNew);CHKERRQ(ierr); 4547bffde78SMichael Lange ierr = PetscMalloc1(numLeaves + numLocalNew, &remotePointsNew);CHKERRQ(ierr); 4557bffde78SMichael Lange for (p = 0; p < numLeaves; ++p) { 4567bffde78SMichael Lange localPointsNew[p] = localPoints[p]; 4577bffde78SMichael Lange remotePointsNew[p].index = remotePoints[p].index; 4587bffde78SMichael Lange remotePointsNew[p].rank = remotePoints[p].rank; 4597bffde78SMichael Lange } 460*f3190fc2SToby Isaac p = numLeaves; 461*f3190fc2SToby Isaac ierr = PetscHashIGetKeys(claimshash, &p, localPointsNew);CHKERRQ(ierr); 462*f3190fc2SToby Isaac ierr = PetscSortInt(numLocalNew,&localPointsNew[numLeaves]);CHKERRQ(ierr); 4637bffde78SMichael Lange for (p = numLeaves; p < numLeaves + numLocalNew; ++p) { 4647bffde78SMichael Lange PetscHashIMap(claimshash, localPointsNew[p], offset); 4657bffde78SMichael Lange remotePointsNew[p] = claims[offset]; 4667bffde78SMichael Lange } 4677bffde78SMichael Lange ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), &sfPointNew);CHKERRQ(ierr); 4687bffde78SMichael Lange ierr = PetscSFSetGraph(sfPointNew, pEnd-pStart, numLeaves+numLocalNew, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr); 4697bffde78SMichael Lange ierr = DMSetPointSF(dm, sfPointNew);CHKERRQ(ierr); 4707bffde78SMichael Lange ierr = PetscSFDestroy(&sfPointNew);CHKERRQ(ierr); 4717bffde78SMichael Lange PetscHashIDestroy(claimshash); 4727bffde78SMichael Lange } 4737bffde78SMichael Lange PetscHashIDestroy(leafhash); 4747bffde78SMichael Lange ierr = PetscHashIJDestroy(&roothash);CHKERRQ(ierr); 4757bffde78SMichael Lange ierr = PetscSectionDestroy(&candidateSection);CHKERRQ(ierr); 4767bffde78SMichael Lange ierr = PetscSectionDestroy(&candidateSectionRemote);CHKERRQ(ierr); 4777bffde78SMichael Lange ierr = PetscSectionDestroy(&claimSection);CHKERRQ(ierr); 4787bffde78SMichael Lange ierr = PetscFree(candidates);CHKERRQ(ierr); 4797bffde78SMichael Lange ierr = PetscFree(candidatesRemote);CHKERRQ(ierr); 4807bffde78SMichael Lange ierr = PetscFree(claims);CHKERRQ(ierr); 48125afeb17SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr); 4827bffde78SMichael Lange PetscFunctionReturn(0); 4837bffde78SMichael Lange } 4847bffde78SMichael Lange 4850c795ddcSMatthew G. Knepley /*@C 48680330477SMatthew G. Knepley DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc. 48780330477SMatthew G. Knepley 48880330477SMatthew G. Knepley Collective on DM 48980330477SMatthew G. Knepley 490e56d480eSMatthew G. Knepley Input Parameters: 491e56d480eSMatthew G. Knepley + dm - The DMPlex object with only cells and vertices 492e56d480eSMatthew G. Knepley - dmInt - If NULL a new DM is created, otherwise the interpolated DM is put into the given DM 49380330477SMatthew G. Knepley 49480330477SMatthew G. Knepley Output Parameter: 4954e3744c5SMatthew G. Knepley . dmInt - The complete DMPlex object 49680330477SMatthew G. Knepley 49780330477SMatthew G. Knepley Level: intermediate 49880330477SMatthew G. Knepley 49980330477SMatthew G. Knepley .keywords: mesh 5004e3744c5SMatthew G. Knepley .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList() 50180330477SMatthew G. Knepley @*/ 5029f074e33SMatthew G Knepley PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt) 5039f074e33SMatthew G Knepley { 5049a5b13a2SMatthew G Knepley DM idm, odm = dm; 5057bffde78SMichael Lange PetscSF sfPoint; 5062e1b13c2SMatthew G. Knepley PetscInt depth, dim, d; 5079f074e33SMatthew G Knepley PetscErrorCode ierr; 5089f074e33SMatthew G Knepley 5099f074e33SMatthew G Knepley PetscFunctionBegin; 510a72f3261SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr); 5112e1b13c2SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 512c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 51376b791aaSMatthew G. Knepley if (dim <= 1) { 51476b791aaSMatthew G. Knepley ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 51576b791aaSMatthew G. Knepley idm = dm; 51676b791aaSMatthew G. Knepley } 5179a5b13a2SMatthew G Knepley for (d = 1; d < dim; ++d) { 5189a5b13a2SMatthew G Knepley /* Create interpolated mesh */ 519e56d480eSMatthew G. Knepley if ((d == dim-1) && *dmInt) {idm = *dmInt;} 520e56d480eSMatthew G. Knepley else {ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr);} 5219a5b13a2SMatthew G Knepley ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr); 522c73cfb54SMatthew G. Knepley ierr = DMSetDimension(idm, dim);CHKERRQ(ierr); 5237bffde78SMichael Lange if (depth > 0) { 5247bffde78SMichael Lange ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr); 5257bffde78SMichael Lange ierr = DMGetPointSF(odm, &sfPoint);CHKERRQ(ierr); 5267bffde78SMichael Lange ierr = DMPlexInterpolatePointSF(idm, sfPoint, depth);CHKERRQ(ierr); 5277bffde78SMichael Lange } 5289a5b13a2SMatthew G Knepley if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);} 5299a5b13a2SMatthew G Knepley odm = idm; 5309f074e33SMatthew G Knepley } 5319a5b13a2SMatthew G Knepley *dmInt = idm; 532a72f3261SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr); 5339f074e33SMatthew G Knepley PetscFunctionReturn(0); 5349f074e33SMatthew G Knepley } 53507b0a1fcSMatthew G Knepley 53680330477SMatthew G. Knepley /*@ 53780330477SMatthew G. Knepley DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices 53880330477SMatthew G. Knepley 53980330477SMatthew G. Knepley Collective on DM 54080330477SMatthew G. Knepley 54180330477SMatthew G. Knepley Input Parameter: 54280330477SMatthew G. Knepley . dmA - The DMPlex object with initial coordinates 54380330477SMatthew G. Knepley 54480330477SMatthew G. Knepley Output Parameter: 54580330477SMatthew G. Knepley . dmB - The DMPlex object with copied coordinates 54680330477SMatthew G. Knepley 54780330477SMatthew G. Knepley Level: intermediate 54880330477SMatthew G. Knepley 54980330477SMatthew G. Knepley Note: This is typically used when adding pieces other than vertices to a mesh 55080330477SMatthew G. Knepley 55180330477SMatthew G. Knepley .keywords: mesh 55265f90189SJed Brown .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMGetCoordinateSection() 55380330477SMatthew G. Knepley @*/ 55407b0a1fcSMatthew G Knepley PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB) 55507b0a1fcSMatthew G Knepley { 55607b0a1fcSMatthew G Knepley Vec coordinatesA, coordinatesB; 55707b0a1fcSMatthew G Knepley PetscSection coordSectionA, coordSectionB; 55807b0a1fcSMatthew G Knepley PetscScalar *coordsA, *coordsB; 55907b0a1fcSMatthew G Knepley PetscInt spaceDim, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d; 56007b0a1fcSMatthew G Knepley PetscErrorCode ierr; 56107b0a1fcSMatthew G Knepley 56207b0a1fcSMatthew G Knepley PetscFunctionBegin; 56376b791aaSMatthew G. Knepley if (dmA == dmB) PetscFunctionReturn(0); 56407b0a1fcSMatthew G Knepley ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr); 56507b0a1fcSMatthew G Knepley ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr); 56607b0a1fcSMatthew 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); 56769d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr); 56869d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr); 569972bc18aSToby Isaac if (coordSectionA == coordSectionB) PetscFunctionReturn(0); 570df26b574SMatthew G. Knepley if (!coordSectionB) { 571df26b574SMatthew G. Knepley PetscInt dim; 572df26b574SMatthew G. Knepley 573df26b574SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) coordSectionA), &coordSectionB);CHKERRQ(ierr); 574df26b574SMatthew G. Knepley ierr = DMGetCoordinateDim(dmA, &dim);CHKERRQ(ierr); 575df26b574SMatthew G. Knepley ierr = DMSetCoordinateSection(dmB, dim, coordSectionB);CHKERRQ(ierr); 576df26b574SMatthew G. Knepley ierr = PetscObjectDereference((PetscObject) coordSectionB);CHKERRQ(ierr); 577df26b574SMatthew G. Knepley } 57807b0a1fcSMatthew G Knepley ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr); 57907b0a1fcSMatthew G Knepley ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr); 58007b0a1fcSMatthew G Knepley ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr); 58107b0a1fcSMatthew G Knepley ierr = PetscSectionSetChart(coordSectionB, vStartB, vEndB);CHKERRQ(ierr); 58207b0a1fcSMatthew G Knepley for (v = vStartB; v < vEndB; ++v) { 58307b0a1fcSMatthew G Knepley ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr); 58407b0a1fcSMatthew G Knepley ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr); 58507b0a1fcSMatthew G Knepley } 58607b0a1fcSMatthew G Knepley ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr); 58707b0a1fcSMatthew G Knepley ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr); 58807b0a1fcSMatthew G Knepley ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr); 5898b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinatesB);CHKERRQ(ierr); 59007b0a1fcSMatthew G Knepley ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr); 59107b0a1fcSMatthew G Knepley ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr); 5922eb5907fSJed Brown ierr = VecSetType(coordinatesB,VECSTANDARD);CHKERRQ(ierr); 59307b0a1fcSMatthew G Knepley ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr); 59407b0a1fcSMatthew G Knepley ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr); 59507b0a1fcSMatthew G Knepley for (v = 0; v < vEndB-vStartB; ++v) { 59607b0a1fcSMatthew G Knepley for (d = 0; d < spaceDim; ++d) { 59707b0a1fcSMatthew G Knepley coordsB[v*spaceDim+d] = coordsA[v*spaceDim+d]; 59807b0a1fcSMatthew G Knepley } 59907b0a1fcSMatthew G Knepley } 60007b0a1fcSMatthew G Knepley ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr); 60107b0a1fcSMatthew G Knepley ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr); 60207b0a1fcSMatthew G Knepley ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr); 60307b0a1fcSMatthew G Knepley ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr); 60407b0a1fcSMatthew G Knepley PetscFunctionReturn(0); 60507b0a1fcSMatthew G Knepley } 6065c386225SMatthew G. Knepley 6074e3744c5SMatthew G. Knepley /*@ 6084e3744c5SMatthew G. Knepley DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh 6094e3744c5SMatthew G. Knepley 6104e3744c5SMatthew G. Knepley Collective on DM 6114e3744c5SMatthew G. Knepley 6124e3744c5SMatthew G. Knepley Input Parameter: 6134e3744c5SMatthew G. Knepley . dm - The complete DMPlex object 6144e3744c5SMatthew G. Knepley 6154e3744c5SMatthew G. Knepley Output Parameter: 6164e3744c5SMatthew G. Knepley . dmUnint - The DMPlex object with only cells and vertices 6174e3744c5SMatthew G. Knepley 6184e3744c5SMatthew G. Knepley Level: intermediate 6194e3744c5SMatthew G. Knepley 6204e3744c5SMatthew G. Knepley .keywords: mesh 6214e3744c5SMatthew G. Knepley .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList() 6224e3744c5SMatthew G. Knepley @*/ 6234e3744c5SMatthew G. Knepley PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint) 6244e3744c5SMatthew G. Knepley { 6254e3744c5SMatthew G. Knepley DM udm; 626595d4abbSMatthew G. Knepley PetscInt dim, vStart, vEnd, cStart, cEnd, c, maxConeSize = 0, *cone; 6274e3744c5SMatthew G. Knepley PetscErrorCode ierr; 6284e3744c5SMatthew G. Knepley 6294e3744c5SMatthew G. Knepley PetscFunctionBegin; 630c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 6314e3744c5SMatthew G. Knepley if (dim <= 1) { 6324e3744c5SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 633595d4abbSMatthew G. Knepley *dmUnint = dm; 634595d4abbSMatthew G. Knepley PetscFunctionReturn(0); 6354e3744c5SMatthew G. Knepley } 636595d4abbSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 637595d4abbSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 6384e3744c5SMatthew G. Knepley ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr); 6394e3744c5SMatthew G. Knepley ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr); 640c73cfb54SMatthew G. Knepley ierr = DMSetDimension(udm, dim);CHKERRQ(ierr); 6414e3744c5SMatthew G. Knepley ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr); 6424e3744c5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 643595d4abbSMatthew G. Knepley PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 6444e3744c5SMatthew G. Knepley 6454e3744c5SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 6464e3744c5SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 6474e3744c5SMatthew G. Knepley const PetscInt p = closure[cl]; 6484e3744c5SMatthew G. Knepley 6494e3744c5SMatthew G. Knepley if ((p >= vStart) && (p < vEnd)) ++coneSize; 6504e3744c5SMatthew G. Knepley } 6514e3744c5SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 6524e3744c5SMatthew G. Knepley ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr); 653595d4abbSMatthew G. Knepley maxConeSize = PetscMax(maxConeSize, coneSize); 6544e3744c5SMatthew G. Knepley } 6554e3744c5SMatthew G. Knepley ierr = DMSetUp(udm);CHKERRQ(ierr); 656785e854fSJed Brown ierr = PetscMalloc1(maxConeSize, &cone);CHKERRQ(ierr); 6574e3744c5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 658595d4abbSMatthew G. Knepley PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 6594e3744c5SMatthew G. Knepley 6604e3744c5SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 6614e3744c5SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 6624e3744c5SMatthew G. Knepley const PetscInt p = closure[cl]; 6634e3744c5SMatthew G. Knepley 6644e3744c5SMatthew G. Knepley if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p; 6654e3744c5SMatthew G. Knepley } 6664e3744c5SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 6674e3744c5SMatthew G. Knepley ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr); 6684e3744c5SMatthew G. Knepley } 6694e3744c5SMatthew G. Knepley ierr = PetscFree(cone);CHKERRQ(ierr); 6704e3744c5SMatthew G. Knepley ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr); 6714e3744c5SMatthew G. Knepley ierr = DMPlexStratify(udm);CHKERRQ(ierr); 6725c7de58aSMatthew G. Knepley /* Reduce SF */ 6735c7de58aSMatthew G. Knepley { 6745c7de58aSMatthew G. Knepley PetscSF sfPoint, sfPointUn; 6755c7de58aSMatthew G. Knepley const PetscSFNode *remotePoints; 6765c7de58aSMatthew G. Knepley const PetscInt *localPoints; 6775c7de58aSMatthew G. Knepley PetscSFNode *remotePointsUn; 6785c7de58aSMatthew G. Knepley PetscInt *localPointsUn; 6795c7de58aSMatthew G. Knepley PetscInt vEnd, numRoots, numLeaves, l; 6805c7de58aSMatthew G. Knepley PetscInt numLeavesUn = 0, n = 0; 6815c7de58aSMatthew G. Knepley PetscErrorCode ierr; 6825c7de58aSMatthew G. Knepley 6835c7de58aSMatthew G. Knepley /* Get original SF information */ 6845c7de58aSMatthew G. Knepley ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 6855c7de58aSMatthew G. Knepley ierr = DMGetPointSF(udm, &sfPointUn);CHKERRQ(ierr); 6865c7de58aSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, NULL, &vEnd);CHKERRQ(ierr); 6875c7de58aSMatthew G. Knepley ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 6885c7de58aSMatthew G. Knepley /* Allocate space for cells and vertices */ 6895c7de58aSMatthew G. Knepley for (l = 0; l < numLeaves; ++l) if (localPoints[l] < vEnd) numLeavesUn++; 6905c7de58aSMatthew G. Knepley /* Fill in leaves */ 6915c7de58aSMatthew G. Knepley if (vEnd >= 0) { 6925c7de58aSMatthew G. Knepley ierr = PetscMalloc1(numLeavesUn, &remotePointsUn);CHKERRQ(ierr); 6935c7de58aSMatthew G. Knepley ierr = PetscMalloc1(numLeavesUn, &localPointsUn);CHKERRQ(ierr); 6945c7de58aSMatthew G. Knepley for (l = 0; l < numLeaves; l++) { 6955c7de58aSMatthew G. Knepley if (localPoints[l] < vEnd) { 6965c7de58aSMatthew G. Knepley localPointsUn[n] = localPoints[l]; 6975c7de58aSMatthew G. Knepley remotePointsUn[n].rank = remotePoints[l].rank; 6985c7de58aSMatthew G. Knepley remotePointsUn[n].index = remotePoints[l].index; 6995c7de58aSMatthew G. Knepley ++n; 7005c7de58aSMatthew G. Knepley } 7015c7de58aSMatthew G. Knepley } 7025c7de58aSMatthew G. Knepley if (n != numLeavesUn) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent number of leaves %d != %d", n, numLeavesUn); 7035c7de58aSMatthew G. Knepley ierr = PetscSFSetGraph(sfPointUn, vEnd, numLeavesUn, localPointsUn, PETSC_OWN_POINTER, remotePointsUn, PETSC_OWN_POINTER);CHKERRQ(ierr); 7045c7de58aSMatthew G. Knepley } 7055c7de58aSMatthew G. Knepley } 7064e3744c5SMatthew G. Knepley *dmUnint = udm; 7074e3744c5SMatthew G. Knepley PetscFunctionReturn(0); 7084e3744c5SMatthew G. Knepley } 709