xref: /petsc/src/dm/impls/plex/plexinterpolate.c (revision 7bffde78dbb5e84675fb811e3c60192f07394e71)
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 #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) {
54c49d9212SMatthew G. Knepley   case 1:
55c49d9212SMatthew G. Knepley     switch (coneSize) {
56c49d9212SMatthew G. Knepley     case 2:
57c49d9212SMatthew G. Knepley       if (faces) {
58c49d9212SMatthew G. Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
59c49d9212SMatthew G. Knepley         *faces = facesTmp;
60c49d9212SMatthew G. Knepley       }
61c49d9212SMatthew G. Knepley       if (numFaces) *numFaces         = 2;
62c49d9212SMatthew G. Knepley       if (faceSize) *faceSize         = 1;
63c49d9212SMatthew G. Knepley       break;
64c49d9212SMatthew G. Knepley     default:
65c49d9212SMatthew G. Knepley       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
66c49d9212SMatthew G. Knepley     }
67c49d9212SMatthew 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;
155c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(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;
163dcca6d9dSJed Brown   ierr = PetscMalloc2(depth+1,&pStart,depth+1,&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   for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) {
1769f074e33SMatthew G Knepley     const PetscInt *cellFaces;
177735a0255SMatthew G. Knepley     PetscInt        numCellFaces, faceSize, cf;
1789f074e33SMatthew G Knepley 
1799a5b13a2SMatthew G Knepley     ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
1809a5b13a2SMatthew G Knepley     if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll);
1819f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
1829a5b13a2SMatthew G Knepley       const PetscInt   *cellFace = &cellFaces[cf*faceSize];
1839a5b13a2SMatthew G Knepley       PetscHashIJKLKey  key;
184735a0255SMatthew G. Knepley       PetscHashIJKLIter missing, iter;
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;
193302440fdSBarry Smith         ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr);
1949f074e33SMatthew G Knepley       }
195735a0255SMatthew G. Knepley       ierr = PetscHashIJKLPut(faceTable, key, &missing, &iter);CHKERRQ(ierr);
196735a0255SMatthew G. Knepley       if (missing) {ierr = PetscHashIJKLSet(faceTable, iter, face++);CHKERRQ(ierr);}
1979a5b13a2SMatthew G Knepley     }
198439ece16SMatthew G. Knepley     ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
1999a5b13a2SMatthew G Knepley   }
2009a5b13a2SMatthew G Knepley   pEnd[faceDepth] = face;
2019a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr);
2029a5b13a2SMatthew G Knepley   /* Count new points */
2039a5b13a2SMatthew G Knepley   for (d = 0; d <= depth; ++d) {
2049a5b13a2SMatthew G Knepley     numPoints += pEnd[d]-pStart[d];
2059a5b13a2SMatthew G Knepley   }
2069a5b13a2SMatthew G Knepley   ierr = DMPlexSetChart(idm, 0, numPoints);CHKERRQ(ierr);
2079a5b13a2SMatthew G Knepley   /* Set cone sizes */
2089a5b13a2SMatthew G Knepley   for (d = 0; d <= depth; ++d) {
2099a5b13a2SMatthew G Knepley     PetscInt coneSize, p;
2109f074e33SMatthew G Knepley 
2119a5b13a2SMatthew G Knepley     if (d == faceDepth) {
2129a5b13a2SMatthew G Knepley       for (p = pStart[d]; p < pEnd[d]; ++p) {
2139a5b13a2SMatthew G Knepley         /* I see no way to do this if we admit faces of different shapes */
2149a5b13a2SMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr);
2159a5b13a2SMatthew G Knepley       }
216a014044eSMatthew G Knepley     } else if (d == cellDepth) {
217a014044eSMatthew G Knepley       for (p = pStart[d]; p < pEnd[d]; ++p) {
218a014044eSMatthew G Knepley         /* Number of cell faces may be different from number of cell vertices*/
219a014044eSMatthew G Knepley         ierr = DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);CHKERRQ(ierr);
220a014044eSMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr);
221a014044eSMatthew G Knepley       }
2229a5b13a2SMatthew G Knepley     } else {
2239a5b13a2SMatthew G Knepley       for (p = pStart[d]; p < pEnd[d]; ++p) {
2249a5b13a2SMatthew G Knepley         ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
2259a5b13a2SMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr);
2269f074e33SMatthew G Knepley       }
2279f074e33SMatthew G Knepley     }
2289f074e33SMatthew G Knepley   }
2299f074e33SMatthew G Knepley   ierr = DMSetUp(idm);CHKERRQ(ierr);
2309f074e33SMatthew G Knepley   /* Get face cones from subsets of cell vertices */
2319a5b13a2SMatthew 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);
2329a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr);
2339a5b13a2SMatthew G Knepley   for (d = depth; d > cellDepth; --d) {
2349a5b13a2SMatthew G Knepley     const PetscInt *cone;
2359a5b13a2SMatthew G Knepley     PetscInt        p;
2369a5b13a2SMatthew G Knepley 
2379a5b13a2SMatthew G Knepley     for (p = pStart[d]; p < pEnd[d]; ++p) {
2389a5b13a2SMatthew G Knepley       ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
2399a5b13a2SMatthew G Knepley       ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr);
2409a5b13a2SMatthew G Knepley       ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr);
2419a5b13a2SMatthew G Knepley       ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr);
2429f074e33SMatthew G Knepley     }
2439a5b13a2SMatthew G Knepley   }
2449a5b13a2SMatthew G Knepley   for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) {
2459f074e33SMatthew G Knepley     const PetscInt *cellFaces;
246735a0255SMatthew G. Knepley     PetscInt        numCellFaces, faceSize, cf;
2479f074e33SMatthew G Knepley 
2489a5b13a2SMatthew G Knepley     ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
2499a5b13a2SMatthew G Knepley     if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll);
2509f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
2519a5b13a2SMatthew G Knepley       const PetscInt  *cellFace = &cellFaces[cf*faceSize];
2529a5b13a2SMatthew G Knepley       PetscHashIJKLKey key;
253735a0255SMatthew G. Knepley       PetscHashIJKLIter missing, iter;
2549f074e33SMatthew G Knepley 
2559a5b13a2SMatthew G Knepley       if (faceSize == 2) {
2569a5b13a2SMatthew G Knepley         key.i = PetscMin(cellFace[0], cellFace[1]);
2579a5b13a2SMatthew G Knepley         key.j = PetscMax(cellFace[0], cellFace[1]);
2587fcd4ef0SJed Brown         key.k = 0;
2597fcd4ef0SJed Brown         key.l = 0;
2609a5b13a2SMatthew G Knepley       } else {
2619a5b13a2SMatthew G Knepley         key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0;
262302440fdSBarry Smith         ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr);
2639f074e33SMatthew G Knepley       }
264735a0255SMatthew G. Knepley       ierr = PetscHashIJKLPut(faceTable, key, &missing, &iter);CHKERRQ(ierr);
265735a0255SMatthew G. Knepley       if (missing) {
2669a5b13a2SMatthew G Knepley         ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr);
267735a0255SMatthew G. Knepley         ierr = PetscHashIJKLSet(faceTable, iter, face);CHKERRQ(ierr);
268735a0255SMatthew G. Knepley         ierr = DMPlexInsertCone(idm, c, cf, face++);CHKERRQ(ierr);
2699a5b13a2SMatthew G Knepley       } else {
2709a5b13a2SMatthew G Knepley         const PetscInt *cone;
271735a0255SMatthew G. Knepley         PetscInt        coneSize, ornt, i, j, f;
2729f074e33SMatthew G Knepley 
273735a0255SMatthew G. Knepley         ierr = PetscHashIJKLGet(faceTable, iter, &f);CHKERRQ(ierr);
2749a5b13a2SMatthew G Knepley         ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr);
2752e1b13c2SMatthew G. Knepley         /* Orient face: Do not allow reverse orientation at the first vertex */
2769f074e33SMatthew G Knepley         ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr);
2779f074e33SMatthew G Knepley         ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr);
2789a5b13a2SMatthew 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);
2799a5b13a2SMatthew G Knepley         /* - First find the initial vertex */
2809a5b13a2SMatthew G Knepley         for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break;
2819a5b13a2SMatthew G Knepley         /* - Try forward comparison */
2829a5b13a2SMatthew G Knepley         for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break;
2839a5b13a2SMatthew G Knepley         if (j == faceSize) {
2849a5b13a2SMatthew G Knepley           if ((faceSize == 2) && (i == 1)) ornt = -2;
2859a5b13a2SMatthew G Knepley           else                             ornt = i;
2869a5b13a2SMatthew G Knepley         } else {
2879a5b13a2SMatthew G Knepley           /* - Try backward comparison */
2889a5b13a2SMatthew G Knepley           for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break;
2892e1b13c2SMatthew G. Knepley           if (j == faceSize) {
2902e1b13c2SMatthew G. Knepley             if (i == 0) ornt = -faceSize;
291dc1a705cSMatthew G. Knepley             else        ornt = -i;
2922e1b13c2SMatthew G. Knepley           } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation");
2939a5b13a2SMatthew G Knepley         }
2949a5b13a2SMatthew G Knepley         ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr);
2959f074e33SMatthew G Knepley       }
2969f074e33SMatthew G Knepley     }
297439ece16SMatthew G. Knepley     ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
2989f074e33SMatthew G Knepley   }
2999a5b13a2SMatthew 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]);
300c907b753SJed Brown   ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr);
3019a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr);
3026551a8c7SMatthew G. Knepley   ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr);
3039a5b13a2SMatthew G Knepley   ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr);
3049a5b13a2SMatthew G Knepley   ierr = DMPlexStratify(idm);CHKERRQ(ierr);
3059f074e33SMatthew G Knepley   PetscFunctionReturn(0);
3069f074e33SMatthew G Knepley }
3079f074e33SMatthew G Knepley 
3089f074e33SMatthew G Knepley #undef __FUNCT__
309*7bffde78SMichael Lange #define __FUNCT__ "DMPlexInterpolatePointSF"
310*7bffde78SMichael Lange /* This interpolates the PointSF in parallel following local interpolation */
311*7bffde78SMichael Lange static PetscErrorCode DMPlexInterpolatePointSF(DM dm, PetscSF pointSF, PetscInt depth)
312*7bffde78SMichael Lange {
313*7bffde78SMichael Lange   PetscMPIInt        rank;
314*7bffde78SMichael Lange   PetscInt           p, c, d, dof, offset;
315*7bffde78SMichael Lange   PetscInt           numLeaves, numRoots, candidatesSize, candidatesRemoteSize, numLeavesInverse;
316*7bffde78SMichael Lange   const PetscInt    *localPoints;
317*7bffde78SMichael Lange   const PetscSFNode *remotePoints;
318*7bffde78SMichael Lange   PetscSFNode       *candidates, *candidatesRemote, *claims;
319*7bffde78SMichael Lange   PetscSection       candidateSection, candidateSectionRemote, claimSection;
320*7bffde78SMichael Lange   PetscHashI         leafhash;
321*7bffde78SMichael Lange   PetscHashIJ        roothash;
322*7bffde78SMichael Lange   PetscHashIJKey     key;
323*7bffde78SMichael Lange   PetscErrorCode     ierr;
324*7bffde78SMichael Lange 
325*7bffde78SMichael Lange   PetscFunctionBegin;
326*7bffde78SMichael Lange   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
327*7bffde78SMichael Lange   ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr);
328*7bffde78SMichael Lange   /* Build hashes of points in the SF for efficient lookup */
329*7bffde78SMichael Lange   PetscHashICreate(leafhash);
330*7bffde78SMichael Lange   PetscHashIJCreate(&roothash);
331*7bffde78SMichael Lange   ierr = PetscHashIJSetMultivalued(roothash, PETSC_FALSE);CHKERRQ(ierr);
332*7bffde78SMichael Lange   for (p = 0; p < numLeaves; ++p) {
333*7bffde78SMichael Lange     PetscHashIAdd(leafhash, localPoints[p], p);
334*7bffde78SMichael Lange     key.i = remotePoints[p].index; key.j = remotePoints[p].rank;
335*7bffde78SMichael Lange     PetscHashIJAdd(roothash, key, p);
336*7bffde78SMichael Lange   }
337*7bffde78SMichael Lange   /* Build a section / SFNode array of candidate points in the single-level adjacency of leaves,
338*7bffde78SMichael Lange      where each candidate is defined by the root entry for the other vertex that defines the edge. */
339*7bffde78SMichael Lange   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSection);CHKERRQ(ierr);
340*7bffde78SMichael Lange   ierr = PetscSectionSetChart(candidateSection, 0, numRoots);CHKERRQ(ierr);
341*7bffde78SMichael Lange   {
342*7bffde78SMichael Lange     PetscInt leaf, root, idx, a, adjSize, *adj = NULL;
343*7bffde78SMichael Lange     for (p = 0; p < numLeaves; ++p) {
344*7bffde78SMichael Lange       PetscInt adjSize = PETSC_DETERMINE;
345*7bffde78SMichael Lange       ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr);
346*7bffde78SMichael Lange       for (a = 0; a < adjSize; ++a) {
347*7bffde78SMichael Lange         PetscHashIMap(leafhash, adj[a], leaf);
348*7bffde78SMichael Lange         if (leaf >= 0) {ierr = PetscSectionAddDof(candidateSection, localPoints[p], 1);CHKERRQ(ierr);}
349*7bffde78SMichael Lange       }
350*7bffde78SMichael Lange     }
351*7bffde78SMichael Lange     ierr = PetscSectionSetUp(candidateSection);CHKERRQ(ierr);
352*7bffde78SMichael Lange     ierr = PetscSectionGetStorageSize(candidateSection, &candidatesSize);CHKERRQ(ierr);
353*7bffde78SMichael Lange     ierr = PetscMalloc1(candidatesSize, &candidates);CHKERRQ(ierr);
354*7bffde78SMichael Lange     for (p = 0; p < numLeaves; ++p) {
355*7bffde78SMichael Lange       PetscInt adjSize = PETSC_DETERMINE;
356*7bffde78SMichael Lange       ierr = PetscSectionGetOffset(candidateSection, localPoints[p], &offset);CHKERRQ(ierr);
357*7bffde78SMichael Lange       ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr);
358*7bffde78SMichael Lange       for (idx = 0, a = 0; a < adjSize; ++a) {
359*7bffde78SMichael Lange         PetscHashIMap(leafhash, adj[a], root);
360*7bffde78SMichael Lange         if (root >= 0) candidates[offset+idx++] = remotePoints[root];
361*7bffde78SMichael Lange       }
362*7bffde78SMichael Lange     }
363*7bffde78SMichael Lange     ierr = PetscFree(adj);CHKERRQ(ierr);
364*7bffde78SMichael Lange   }
365*7bffde78SMichael Lange   /* Gather candidate section / array pair into the root partition via inverse(multi(pointSF)). */
366*7bffde78SMichael Lange   {
367*7bffde78SMichael Lange     PetscSF   sfMulti, sfInverse, sfCandidates;
368*7bffde78SMichael Lange     PetscInt *remoteOffsets;
369*7bffde78SMichael Lange     ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr);
370*7bffde78SMichael Lange     ierr = PetscSFCreateInverseSF(sfMulti, &sfInverse);CHKERRQ(ierr);
371*7bffde78SMichael Lange     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSectionRemote);CHKERRQ(ierr);
372*7bffde78SMichael Lange     ierr = PetscSFDistributeSection(sfInverse, candidateSection, &remoteOffsets, candidateSectionRemote);CHKERRQ(ierr);
373*7bffde78SMichael Lange     ierr = PetscSFCreateSectionSF(sfInverse, candidateSection, remoteOffsets, candidateSectionRemote, &sfCandidates);CHKERRQ(ierr);
374*7bffde78SMichael Lange     ierr = PetscSectionGetStorageSize(candidateSectionRemote, &candidatesRemoteSize);CHKERRQ(ierr);
375*7bffde78SMichael Lange     ierr = PetscMalloc1(candidatesRemoteSize, &candidatesRemote);CHKERRQ(ierr);
376*7bffde78SMichael Lange     ierr = PetscSFBcastBegin(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr);
377*7bffde78SMichael Lange     ierr = PetscSFBcastEnd(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr);
378*7bffde78SMichael Lange     ierr = PetscSFDestroy(&sfInverse);CHKERRQ(ierr);
379*7bffde78SMichael Lange     ierr = PetscSFDestroy(&sfCandidates);CHKERRQ(ierr);
380*7bffde78SMichael Lange     ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
381*7bffde78SMichael Lange   }
382*7bffde78SMichael Lange   /* Walk local roots and check for each remote candidate whether we know all required points,
383*7bffde78SMichael Lange      either from owning it or having a root entry in the point SF. If we do we place a claim
384*7bffde78SMichael Lange      by replacing the vertex number with our edge ID. */
385*7bffde78SMichael Lange   {
386*7bffde78SMichael Lange     PetscInt        idx, root, joinSize, vertices[2];
387*7bffde78SMichael Lange     const PetscInt *rootdegree, *join = NULL;
388*7bffde78SMichael Lange     ierr = PetscSFComputeDegreeBegin(pointSF, &rootdegree);CHKERRQ(ierr);
389*7bffde78SMichael Lange     ierr = PetscSFComputeDegreeEnd(pointSF, &rootdegree);CHKERRQ(ierr);
390*7bffde78SMichael Lange     /* Loop remote edge connections and put in a claim if both vertices are known */
391*7bffde78SMichael Lange     for (idx = 0, p = 0; p < numRoots; ++p) {
392*7bffde78SMichael Lange       for (d = 0; d < rootdegree[p]; ++d) {
393*7bffde78SMichael Lange         ierr = PetscSectionGetDof(candidateSectionRemote, idx, &dof);CHKERRQ(ierr);
394*7bffde78SMichael Lange         ierr = PetscSectionGetOffset(candidateSectionRemote, idx, &offset);CHKERRQ(ierr);
395*7bffde78SMichael Lange         for (c = 0; c < dof; ++c) {
396*7bffde78SMichael Lange           /* We own both vertices, so we claim the edge by replacing vertex with edge */
397*7bffde78SMichael Lange           if (candidatesRemote[offset+c].rank == rank) {
398*7bffde78SMichael Lange             vertices[0] = p; vertices[1] = candidatesRemote[offset+c].index;
399*7bffde78SMichael Lange             ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
400*7bffde78SMichael Lange             if (joinSize == 1) candidatesRemote[offset+c].index = join[0];
401*7bffde78SMichael Lange             ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
402*7bffde78SMichael Lange             continue;
403*7bffde78SMichael Lange           }
404*7bffde78SMichael Lange           /* If we own one vertex and share a root with the other, we claim it */
405*7bffde78SMichael Lange           key.i = candidatesRemote[offset+c].index; key.j = candidatesRemote[offset+c].rank;
406*7bffde78SMichael Lange           PetscHashIJGet(roothash, key, &root);
407*7bffde78SMichael Lange           if (root >= 0) {
408*7bffde78SMichael Lange             vertices[0] = p; vertices[1] = localPoints[root];
409*7bffde78SMichael Lange             ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
410*7bffde78SMichael Lange             if (joinSize == 1) {
411*7bffde78SMichael Lange               candidatesRemote[offset+c].index = join[0];
412*7bffde78SMichael Lange               candidatesRemote[offset+c].rank = rank;
413*7bffde78SMichael Lange             }
414*7bffde78SMichael Lange             ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
415*7bffde78SMichael Lange           }
416*7bffde78SMichael Lange         }
417*7bffde78SMichael Lange         idx++;
418*7bffde78SMichael Lange       }
419*7bffde78SMichael Lange     }
420*7bffde78SMichael Lange   }
421*7bffde78SMichael Lange   /* Push claims back to receiver via the MultiSF and derive new pointSF mapping on receiver */
422*7bffde78SMichael Lange   {
423*7bffde78SMichael Lange     PetscSF         sfMulti, sfClaims, sfPointNew;
424*7bffde78SMichael Lange     PetscHashI      claimshash;
425*7bffde78SMichael Lange     PetscInt        size, pStart, pEnd, root, joinSize, numLocalNew;
426*7bffde78SMichael Lange     PetscInt       *remoteOffsets, *localPointsNew, vertices[2];
427*7bffde78SMichael Lange     const PetscInt *join = NULL;;
428*7bffde78SMichael Lange     PetscSFNode    *remotePointsNew;
429*7bffde78SMichael Lange     ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr);
430*7bffde78SMichael Lange     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &claimSection);CHKERRQ(ierr);
431*7bffde78SMichael Lange     ierr = PetscSFDistributeSection(sfMulti, candidateSectionRemote, &remoteOffsets, claimSection);CHKERRQ(ierr);
432*7bffde78SMichael Lange     ierr = PetscSFCreateSectionSF(sfMulti, candidateSectionRemote, remoteOffsets, claimSection, &sfClaims);CHKERRQ(ierr);
433*7bffde78SMichael Lange     ierr = PetscSectionGetStorageSize(claimSection, &size);CHKERRQ(ierr);
434*7bffde78SMichael Lange     ierr = PetscMalloc1(size, &claims);CHKERRQ(ierr);
435*7bffde78SMichael Lange     ierr = PetscSFBcastBegin(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr);
436*7bffde78SMichael Lange     ierr = PetscSFBcastEnd(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr);
437*7bffde78SMichael Lange     ierr = PetscSFDestroy(&sfClaims);CHKERRQ(ierr);
438*7bffde78SMichael Lange     ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
439*7bffde78SMichael Lange     /* Walk the original section of local supports and add an SF entry for each updated item */
440*7bffde78SMichael Lange     PetscHashICreate(claimshash);
441*7bffde78SMichael Lange     for (p = 0; p < numRoots; ++p) {
442*7bffde78SMichael Lange       ierr = PetscSectionGetDof(candidateSection, p, &dof);CHKERRQ(ierr);
443*7bffde78SMichael Lange       ierr = PetscSectionGetOffset(candidateSection, p, &offset);CHKERRQ(ierr);
444*7bffde78SMichael Lange       for (d = 0; d < dof; ++d) {
445*7bffde78SMichael Lange         if (candidates[offset+d].index != claims[offset+d].index) {
446*7bffde78SMichael Lange           key.i = candidates[offset+d].index; key.j = candidates[offset+d].rank;
447*7bffde78SMichael Lange           PetscHashIJGet(roothash, key, &root);
448*7bffde78SMichael Lange           if (root >= 0) {
449*7bffde78SMichael Lange             vertices[0] = p; vertices[1] = localPoints[root];
450*7bffde78SMichael Lange             ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
451*7bffde78SMichael Lange             if (joinSize == 1) PetscHashIAdd(claimshash, join[0], offset+d);
452*7bffde78SMichael Lange             ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
453*7bffde78SMichael Lange           }
454*7bffde78SMichael Lange         }
455*7bffde78SMichael Lange       }
456*7bffde78SMichael Lange     }
457*7bffde78SMichael Lange     /* Create new pointSF from hashed claims */
458*7bffde78SMichael Lange     PetscHashISize(claimshash, numLocalNew);
459*7bffde78SMichael Lange     ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
460*7bffde78SMichael Lange     ierr = PetscMalloc1(numLeaves + numLocalNew, &localPointsNew);CHKERRQ(ierr);
461*7bffde78SMichael Lange     ierr = PetscMalloc1(numLeaves + numLocalNew, &remotePointsNew);CHKERRQ(ierr);
462*7bffde78SMichael Lange     for (p = 0; p < numLeaves; ++p) {
463*7bffde78SMichael Lange       localPointsNew[p] = localPoints[p];
464*7bffde78SMichael Lange       remotePointsNew[p].index = remotePoints[p].index;
465*7bffde78SMichael Lange       remotePointsNew[p].rank = remotePoints[p].rank;
466*7bffde78SMichael Lange     }
467*7bffde78SMichael Lange     p = numLeaves; ierr = PetscHashIGetKeys(claimshash, &p, localPointsNew);CHKERRQ(ierr);
468*7bffde78SMichael Lange     for (p = numLeaves; p < numLeaves + numLocalNew; ++p) {
469*7bffde78SMichael Lange       PetscHashIMap(claimshash, localPointsNew[p], offset);
470*7bffde78SMichael Lange       remotePointsNew[p] = claims[offset];
471*7bffde78SMichael Lange     }
472*7bffde78SMichael Lange     ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), &sfPointNew);CHKERRQ(ierr);
473*7bffde78SMichael Lange     ierr = PetscSFSetGraph(sfPointNew, pEnd-pStart, numLeaves+numLocalNew, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr);
474*7bffde78SMichael Lange     ierr = DMSetPointSF(dm, sfPointNew);CHKERRQ(ierr);
475*7bffde78SMichael Lange     ierr = PetscSFDestroy(&sfPointNew);CHKERRQ(ierr);
476*7bffde78SMichael Lange     PetscHashIDestroy(claimshash);
477*7bffde78SMichael Lange   }
478*7bffde78SMichael Lange   PetscHashIDestroy(leafhash);
479*7bffde78SMichael Lange   ierr = PetscHashIJDestroy(&roothash);CHKERRQ(ierr);
480*7bffde78SMichael Lange   ierr = PetscSectionDestroy(&candidateSection);CHKERRQ(ierr);
481*7bffde78SMichael Lange   ierr = PetscSectionDestroy(&candidateSectionRemote);CHKERRQ(ierr);
482*7bffde78SMichael Lange   ierr = PetscSectionDestroy(&claimSection);CHKERRQ(ierr);
483*7bffde78SMichael Lange   ierr = PetscFree(candidates);CHKERRQ(ierr);
484*7bffde78SMichael Lange   ierr = PetscFree(candidatesRemote);CHKERRQ(ierr);
485*7bffde78SMichael Lange   ierr = PetscFree(claims);CHKERRQ(ierr);
486*7bffde78SMichael Lange   PetscFunctionReturn(0);
487*7bffde78SMichael Lange }
488*7bffde78SMichael Lange 
489*7bffde78SMichael Lange #undef __FUNCT__
4909f074e33SMatthew G Knepley #define __FUNCT__ "DMPlexInterpolate"
4910c795ddcSMatthew G. Knepley /*@C
49280330477SMatthew G. Knepley   DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc.
49380330477SMatthew G. Knepley 
49480330477SMatthew G. Knepley   Collective on DM
49580330477SMatthew G. Knepley 
496e56d480eSMatthew G. Knepley   Input Parameters:
497e56d480eSMatthew G. Knepley + dm - The DMPlex object with only cells and vertices
498e56d480eSMatthew G. Knepley - dmInt - If NULL a new DM is created, otherwise the interpolated DM is put into the given DM
49980330477SMatthew G. Knepley 
50080330477SMatthew G. Knepley   Output Parameter:
5014e3744c5SMatthew G. Knepley . dmInt - The complete DMPlex object
50280330477SMatthew G. Knepley 
50380330477SMatthew G. Knepley   Level: intermediate
50480330477SMatthew G. Knepley 
50580330477SMatthew G. Knepley .keywords: mesh
5064e3744c5SMatthew G. Knepley .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList()
50780330477SMatthew G. Knepley @*/
5089f074e33SMatthew G Knepley PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt)
5099f074e33SMatthew G Knepley {
5109a5b13a2SMatthew G Knepley   DM             idm, odm = dm;
511*7bffde78SMichael Lange   PetscSF        sfPoint;
5122e1b13c2SMatthew G. Knepley   PetscInt       depth, dim, d;
5139f074e33SMatthew G Knepley   PetscErrorCode ierr;
5149f074e33SMatthew G Knepley 
5159f074e33SMatthew G Knepley   PetscFunctionBegin;
516a72f3261SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr);
5172e1b13c2SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
518c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
51976b791aaSMatthew G. Knepley   if (dim <= 1) {
52076b791aaSMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
52176b791aaSMatthew G. Knepley     idm  = dm;
52276b791aaSMatthew G. Knepley   }
5239a5b13a2SMatthew G Knepley   for (d = 1; d < dim; ++d) {
5249a5b13a2SMatthew G Knepley     /* Create interpolated mesh */
525e56d480eSMatthew G. Knepley     if ((d == dim-1) && *dmInt) {idm  = *dmInt;}
526e56d480eSMatthew G. Knepley     else                        {ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr);}
5279a5b13a2SMatthew G Knepley     ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr);
528c73cfb54SMatthew G. Knepley     ierr = DMSetDimension(idm, dim);CHKERRQ(ierr);
529*7bffde78SMichael Lange     if (depth > 0) {
530*7bffde78SMichael Lange       ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr);
531*7bffde78SMichael Lange       ierr = DMGetPointSF(odm, &sfPoint);CHKERRQ(ierr);
532*7bffde78SMichael Lange       ierr = DMPlexInterpolatePointSF(idm, sfPoint, depth);CHKERRQ(ierr);
533*7bffde78SMichael Lange     }
5349a5b13a2SMatthew G Knepley     if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);}
5359a5b13a2SMatthew G Knepley     odm  = idm;
5369f074e33SMatthew G Knepley   }
5379a5b13a2SMatthew G Knepley   *dmInt = idm;
538a72f3261SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr);
5399f074e33SMatthew G Knepley   PetscFunctionReturn(0);
5409f074e33SMatthew G Knepley }
54107b0a1fcSMatthew G Knepley 
54207b0a1fcSMatthew G Knepley #undef __FUNCT__
54307b0a1fcSMatthew G Knepley #define __FUNCT__ "DMPlexCopyCoordinates"
54480330477SMatthew G. Knepley /*@
54580330477SMatthew G. Knepley   DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices
54680330477SMatthew G. Knepley 
54780330477SMatthew G. Knepley   Collective on DM
54880330477SMatthew G. Knepley 
54980330477SMatthew G. Knepley   Input Parameter:
55080330477SMatthew G. Knepley . dmA - The DMPlex object with initial coordinates
55180330477SMatthew G. Knepley 
55280330477SMatthew G. Knepley   Output Parameter:
55380330477SMatthew G. Knepley . dmB - The DMPlex object with copied coordinates
55480330477SMatthew G. Knepley 
55580330477SMatthew G. Knepley   Level: intermediate
55680330477SMatthew G. Knepley 
55780330477SMatthew G. Knepley   Note: This is typically used when adding pieces other than vertices to a mesh
55880330477SMatthew G. Knepley 
55980330477SMatthew G. Knepley .keywords: mesh
56065f90189SJed Brown .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMGetCoordinateSection()
56180330477SMatthew G. Knepley @*/
56207b0a1fcSMatthew G Knepley PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB)
56307b0a1fcSMatthew G Knepley {
56407b0a1fcSMatthew G Knepley   Vec            coordinatesA, coordinatesB;
56507b0a1fcSMatthew G Knepley   PetscSection   coordSectionA, coordSectionB;
56607b0a1fcSMatthew G Knepley   PetscScalar   *coordsA, *coordsB;
56707b0a1fcSMatthew G Knepley   PetscInt       spaceDim, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d;
56807b0a1fcSMatthew G Knepley   PetscErrorCode ierr;
56907b0a1fcSMatthew G Knepley 
57007b0a1fcSMatthew G Knepley   PetscFunctionBegin;
57176b791aaSMatthew G. Knepley   if (dmA == dmB) PetscFunctionReturn(0);
57207b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr);
57307b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr);
57407b0a1fcSMatthew 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);
57569d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr);
57669d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr);
577df26b574SMatthew G. Knepley   if (!coordSectionB) {
578df26b574SMatthew G. Knepley     PetscInt dim;
579df26b574SMatthew G. Knepley 
580df26b574SMatthew G. Knepley     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) coordSectionA), &coordSectionB);CHKERRQ(ierr);
581df26b574SMatthew G. Knepley     ierr = DMGetCoordinateDim(dmA, &dim);CHKERRQ(ierr);
582df26b574SMatthew G. Knepley     ierr = DMSetCoordinateSection(dmB, dim, coordSectionB);CHKERRQ(ierr);
583df26b574SMatthew G. Knepley     ierr = PetscObjectDereference((PetscObject) coordSectionB);CHKERRQ(ierr);
584df26b574SMatthew G. Knepley   }
58507b0a1fcSMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr);
58607b0a1fcSMatthew G Knepley   ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr);
58707b0a1fcSMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr);
58807b0a1fcSMatthew G Knepley   ierr = PetscSectionSetChart(coordSectionB, vStartB, vEndB);CHKERRQ(ierr);
58907b0a1fcSMatthew G Knepley   for (v = vStartB; v < vEndB; ++v) {
59007b0a1fcSMatthew G Knepley     ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr);
59107b0a1fcSMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr);
59207b0a1fcSMatthew G Knepley   }
59307b0a1fcSMatthew G Knepley   ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr);
59407b0a1fcSMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr);
59507b0a1fcSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr);
59607b0a1fcSMatthew G Knepley   ierr = VecCreate(PetscObjectComm((PetscObject) dmB), &coordinatesB);CHKERRQ(ierr);
59707b0a1fcSMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr);
59807b0a1fcSMatthew G Knepley   ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr);
5992eb5907fSJed Brown   ierr = VecSetType(coordinatesB,VECSTANDARD);CHKERRQ(ierr);
60007b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr);
60107b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr);
60207b0a1fcSMatthew G Knepley   for (v = 0; v < vEndB-vStartB; ++v) {
60307b0a1fcSMatthew G Knepley     for (d = 0; d < spaceDim; ++d) {
60407b0a1fcSMatthew G Knepley       coordsB[v*spaceDim+d] = coordsA[v*spaceDim+d];
60507b0a1fcSMatthew G Knepley     }
60607b0a1fcSMatthew G Knepley   }
60707b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr);
60807b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr);
60907b0a1fcSMatthew G Knepley   ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr);
61007b0a1fcSMatthew G Knepley   ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr);
61107b0a1fcSMatthew G Knepley   PetscFunctionReturn(0);
61207b0a1fcSMatthew G Knepley }
6135c386225SMatthew G. Knepley 
6145c386225SMatthew G. Knepley #undef __FUNCT__
6154e3744c5SMatthew G. Knepley #define __FUNCT__ "DMPlexUninterpolate"
6164e3744c5SMatthew G. Knepley /*@
6174e3744c5SMatthew G. Knepley   DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh
6184e3744c5SMatthew G. Knepley 
6194e3744c5SMatthew G. Knepley   Collective on DM
6204e3744c5SMatthew G. Knepley 
6214e3744c5SMatthew G. Knepley   Input Parameter:
6224e3744c5SMatthew G. Knepley . dm - The complete DMPlex object
6234e3744c5SMatthew G. Knepley 
6244e3744c5SMatthew G. Knepley   Output Parameter:
6254e3744c5SMatthew G. Knepley . dmUnint - The DMPlex object with only cells and vertices
6264e3744c5SMatthew G. Knepley 
6274e3744c5SMatthew G. Knepley   Level: intermediate
6284e3744c5SMatthew G. Knepley 
6294e3744c5SMatthew G. Knepley .keywords: mesh
6304e3744c5SMatthew G. Knepley .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList()
6314e3744c5SMatthew G. Knepley @*/
6324e3744c5SMatthew G. Knepley PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint)
6334e3744c5SMatthew G. Knepley {
6344e3744c5SMatthew G. Knepley   DM             udm;
635595d4abbSMatthew G. Knepley   PetscInt       dim, vStart, vEnd, cStart, cEnd, c, maxConeSize = 0, *cone;
6364e3744c5SMatthew G. Knepley   PetscErrorCode ierr;
6374e3744c5SMatthew G. Knepley 
6384e3744c5SMatthew G. Knepley   PetscFunctionBegin;
639c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
6404e3744c5SMatthew G. Knepley   if (dim <= 1) {
6414e3744c5SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
642595d4abbSMatthew G. Knepley     *dmUnint = dm;
643595d4abbSMatthew G. Knepley     PetscFunctionReturn(0);
6444e3744c5SMatthew G. Knepley   }
645595d4abbSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
646595d4abbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
6474e3744c5SMatthew G. Knepley   ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr);
6484e3744c5SMatthew G. Knepley   ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr);
649c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(udm, dim);CHKERRQ(ierr);
6504e3744c5SMatthew G. Knepley   ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr);
6514e3744c5SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
652595d4abbSMatthew G. Knepley     PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
6534e3744c5SMatthew G. Knepley 
6544e3744c5SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
6554e3744c5SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
6564e3744c5SMatthew G. Knepley       const PetscInt p = closure[cl];
6574e3744c5SMatthew G. Knepley 
6584e3744c5SMatthew G. Knepley       if ((p >= vStart) && (p < vEnd)) ++coneSize;
6594e3744c5SMatthew G. Knepley     }
6604e3744c5SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
6614e3744c5SMatthew G. Knepley     ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr);
662595d4abbSMatthew G. Knepley     maxConeSize = PetscMax(maxConeSize, coneSize);
6634e3744c5SMatthew G. Knepley   }
6644e3744c5SMatthew G. Knepley   ierr = DMSetUp(udm);CHKERRQ(ierr);
665785e854fSJed Brown   ierr = PetscMalloc1(maxConeSize, &cone);CHKERRQ(ierr);
6664e3744c5SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
667595d4abbSMatthew G. Knepley     PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
6684e3744c5SMatthew G. Knepley 
6694e3744c5SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
6704e3744c5SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
6714e3744c5SMatthew G. Knepley       const PetscInt p = closure[cl];
6724e3744c5SMatthew G. Knepley 
6734e3744c5SMatthew G. Knepley       if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p;
6744e3744c5SMatthew G. Knepley     }
6754e3744c5SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
6764e3744c5SMatthew G. Knepley     ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr);
6774e3744c5SMatthew G. Knepley   }
6784e3744c5SMatthew G. Knepley   ierr = PetscFree(cone);CHKERRQ(ierr);
6794e3744c5SMatthew G. Knepley   ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr);
6804e3744c5SMatthew G. Knepley   ierr = DMPlexStratify(udm);CHKERRQ(ierr);
6815c7de58aSMatthew G. Knepley   /* Reduce SF */
6825c7de58aSMatthew G. Knepley   {
6835c7de58aSMatthew G. Knepley     PetscSF            sfPoint, sfPointUn;
6845c7de58aSMatthew G. Knepley     const PetscSFNode *remotePoints;
6855c7de58aSMatthew G. Knepley     const PetscInt    *localPoints;
6865c7de58aSMatthew G. Knepley     PetscSFNode       *remotePointsUn;
6875c7de58aSMatthew G. Knepley     PetscInt          *localPointsUn;
6885c7de58aSMatthew G. Knepley     PetscInt           vEnd, numRoots, numLeaves, l;
6895c7de58aSMatthew G. Knepley     PetscInt           numLeavesUn = 0, n = 0;
6905c7de58aSMatthew G. Knepley     PetscErrorCode     ierr;
6915c7de58aSMatthew G. Knepley 
6925c7de58aSMatthew G. Knepley     /* Get original SF information */
6935c7de58aSMatthew G. Knepley     ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
6945c7de58aSMatthew G. Knepley     ierr = DMGetPointSF(udm, &sfPointUn);CHKERRQ(ierr);
6955c7de58aSMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, 0, NULL, &vEnd);CHKERRQ(ierr);
6965c7de58aSMatthew G. Knepley     ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr);
6975c7de58aSMatthew G. Knepley     /* Allocate space for cells and vertices */
6985c7de58aSMatthew G. Knepley     for (l = 0; l < numLeaves; ++l) if (localPoints[l] < vEnd) numLeavesUn++;
6995c7de58aSMatthew G. Knepley     /* Fill in leaves */
7005c7de58aSMatthew G. Knepley     if (vEnd >= 0) {
7015c7de58aSMatthew G. Knepley       ierr = PetscMalloc1(numLeavesUn, &remotePointsUn);CHKERRQ(ierr);
7025c7de58aSMatthew G. Knepley       ierr = PetscMalloc1(numLeavesUn, &localPointsUn);CHKERRQ(ierr);
7035c7de58aSMatthew G. Knepley       for (l = 0; l < numLeaves; l++) {
7045c7de58aSMatthew G. Knepley         if (localPoints[l] < vEnd) {
7055c7de58aSMatthew G. Knepley           localPointsUn[n]        = localPoints[l];
7065c7de58aSMatthew G. Knepley           remotePointsUn[n].rank  = remotePoints[l].rank;
7075c7de58aSMatthew G. Knepley           remotePointsUn[n].index = remotePoints[l].index;
7085c7de58aSMatthew G. Knepley           ++n;
7095c7de58aSMatthew G. Knepley         }
7105c7de58aSMatthew G. Knepley       }
7115c7de58aSMatthew G. Knepley       if (n != numLeavesUn) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent number of leaves %d != %d", n, numLeavesUn);
7125c7de58aSMatthew G. Knepley       ierr = PetscSFSetGraph(sfPointUn, vEnd, numLeavesUn, localPointsUn, PETSC_OWN_POINTER, remotePointsUn, PETSC_OWN_POINTER);CHKERRQ(ierr);
7135c7de58aSMatthew G. Knepley     }
7145c7de58aSMatthew G. Knepley   }
7154e3744c5SMatthew G. Knepley   *dmUnint = udm;
7164e3744c5SMatthew G. Knepley   PetscFunctionReturn(0);
7174e3744c5SMatthew G. Knepley }
718