xref: /petsc/src/dm/impls/plex/plexinterpolate.c (revision 0bedd6be8d0a0fe1e635535d123d6ee49a3cc19e)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
29f074e33SMatthew G Knepley 
39f074e33SMatthew G Knepley /*
49a5b13a2SMatthew G Knepley   DMPlexGetFaces_Internal - Gets groups of vertices that correspond to faces for the given cell
59f074e33SMatthew G Knepley */
6439ece16SMatthew G. Knepley PetscErrorCode DMPlexGetFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[])
79f074e33SMatthew G Knepley {
89f074e33SMatthew G Knepley   const PetscInt *cone = NULL;
99a5b13a2SMatthew G Knepley   PetscInt        maxConeSize, maxSupportSize, coneSize;
109f074e33SMatthew G Knepley   PetscErrorCode  ierr;
119f074e33SMatthew G Knepley 
129f074e33SMatthew G Knepley   PetscFunctionBegin;
139f074e33SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
149a5b13a2SMatthew G Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
159f074e33SMatthew G Knepley   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
169f074e33SMatthew G Knepley   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
17439ece16SMatthew G. Knepley   ierr = DMPlexGetRawFaces_Internal(dm, dim, coneSize, cone, numFaces, faceSize, faces);CHKERRQ(ierr);
18439ece16SMatthew G. Knepley   PetscFunctionReturn(0);
19439ece16SMatthew G. Knepley }
20439ece16SMatthew G. Knepley 
21439ece16SMatthew G. Knepley /*
22439ece16SMatthew G. Knepley   DMPlexRestoreFaces_Internal - Restores the array
23439ece16SMatthew G. Knepley */
24439ece16SMatthew G. Knepley PetscErrorCode DMPlexRestoreFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[])
25439ece16SMatthew G. Knepley {
26439ece16SMatthew G. Knepley   PetscErrorCode  ierr;
27439ece16SMatthew G. Knepley 
28439ece16SMatthew G. Knepley   PetscFunctionBegin;
2969291d52SBarry Smith   ierr = DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faces);CHKERRQ(ierr);
30439ece16SMatthew G. Knepley   PetscFunctionReturn(0);
31439ece16SMatthew G. Knepley }
32439ece16SMatthew G. Knepley 
33439ece16SMatthew G. Knepley /*
34439ece16SMatthew G. Knepley   DMPlexGetRawFaces_Internal - Gets groups of vertices that correspond to faces for the given cone
35439ece16SMatthew G. Knepley */
36439ece16SMatthew G. Knepley PetscErrorCode DMPlexGetRawFaces_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[])
37439ece16SMatthew G. Knepley {
38439ece16SMatthew G. Knepley   PetscInt       *facesTmp;
39439ece16SMatthew G. Knepley   PetscInt        maxConeSize, maxSupportSize;
40439ece16SMatthew G. Knepley   PetscErrorCode  ierr;
41439ece16SMatthew G. Knepley 
42439ece16SMatthew G. Knepley   PetscFunctionBegin;
43439ece16SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
44439ece16SMatthew G. Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
4569291d52SBarry Smith   if (faces) {ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), MPIU_INT, &facesTmp);CHKERRQ(ierr);}
469f074e33SMatthew G Knepley   switch (dim) {
47c49d9212SMatthew G. Knepley   case 1:
48c49d9212SMatthew G. Knepley     switch (coneSize) {
49c49d9212SMatthew G. Knepley     case 2:
50c49d9212SMatthew G. Knepley       if (faces) {
51c49d9212SMatthew G. Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
52c49d9212SMatthew G. Knepley         *faces = facesTmp;
53c49d9212SMatthew G. Knepley       }
54c49d9212SMatthew G. Knepley       if (numFaces) *numFaces         = 2;
55c49d9212SMatthew G. Knepley       if (faceSize) *faceSize         = 1;
56c49d9212SMatthew G. Knepley       break;
57c49d9212SMatthew G. Knepley     default:
58c49d9212SMatthew G. Knepley       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
59c49d9212SMatthew G. Knepley     }
60c49d9212SMatthew G. Knepley     break;
619f074e33SMatthew G Knepley   case 2:
629f074e33SMatthew G Knepley     switch (coneSize) {
639f074e33SMatthew G Knepley     case 3:
649a5b13a2SMatthew G Knepley       if (faces) {
659a5b13a2SMatthew G Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
669a5b13a2SMatthew G Knepley         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
679a5b13a2SMatthew G Knepley         facesTmp[4] = cone[2]; facesTmp[5] = cone[0];
689a5b13a2SMatthew G Knepley         *faces = facesTmp;
699a5b13a2SMatthew G Knepley       }
709a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 3;
719a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 2;
729f074e33SMatthew G Knepley       break;
739f074e33SMatthew G Knepley     case 4:
749a5b13a2SMatthew G Knepley       /* Vertices follow right hand rule */
759a5b13a2SMatthew G Knepley       if (faces) {
769a5b13a2SMatthew G Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
779a5b13a2SMatthew G Knepley         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
789a5b13a2SMatthew G Knepley         facesTmp[4] = cone[2]; facesTmp[5] = cone[3];
799a5b13a2SMatthew G Knepley         facesTmp[6] = cone[3]; facesTmp[7] = cone[0];
809a5b13a2SMatthew G Knepley         *faces = facesTmp;
819a5b13a2SMatthew G Knepley       }
829a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 4;
839a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 2;
849f074e33SMatthew G Knepley       break;
859f074e33SMatthew G Knepley     default:
869f074e33SMatthew G Knepley       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
879f074e33SMatthew G Knepley     }
889f074e33SMatthew G Knepley     break;
899f074e33SMatthew G Knepley   case 3:
909f074e33SMatthew G Knepley     switch (coneSize) {
919f074e33SMatthew G Knepley     case 3:
929a5b13a2SMatthew G Knepley       if (faces) {
939a5b13a2SMatthew G Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
949a5b13a2SMatthew G Knepley         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
959a5b13a2SMatthew G Knepley         facesTmp[4] = cone[2]; facesTmp[5] = cone[0];
969a5b13a2SMatthew G Knepley         *faces = facesTmp;
979a5b13a2SMatthew G Knepley       }
989a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 3;
999a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 2;
1009f074e33SMatthew G Knepley       break;
1019f074e33SMatthew G Knepley     case 4:
1022e1b13c2SMatthew G. Knepley       /* Vertices of first face follow right hand rule and normal points away from last vertex */
1039a5b13a2SMatthew G Knepley       if (faces) {
1042e1b13c2SMatthew G. Knepley         facesTmp[0] = cone[0]; facesTmp[1]  = cone[1]; facesTmp[2]  = cone[2];
1052e1b13c2SMatthew G. Knepley         facesTmp[3] = cone[0]; facesTmp[4]  = cone[3]; facesTmp[5]  = cone[1];
1062e1b13c2SMatthew G. Knepley         facesTmp[6] = cone[0]; facesTmp[7]  = cone[2]; facesTmp[8]  = cone[3];
1072e1b13c2SMatthew G. Knepley         facesTmp[9] = cone[2]; facesTmp[10] = cone[1]; facesTmp[11] = cone[3];
1089a5b13a2SMatthew G Knepley         *faces = facesTmp;
1099a5b13a2SMatthew G Knepley       }
1109a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 4;
1119a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 3;
1129a5b13a2SMatthew G Knepley       break;
1139a5b13a2SMatthew G Knepley     case 8:
1149a5b13a2SMatthew G Knepley       if (faces) {
11551cfd6a4SMatthew G. Knepley         facesTmp[0]  = cone[0]; facesTmp[1]  = cone[1]; facesTmp[2]  = cone[2]; facesTmp[3]  = cone[3]; /* Bottom */
11651cfd6a4SMatthew G. Knepley         facesTmp[4]  = cone[4]; facesTmp[5]  = cone[5]; facesTmp[6]  = cone[6]; facesTmp[7]  = cone[7]; /* Top */
11751cfd6a4SMatthew G. Knepley         facesTmp[8]  = cone[0]; facesTmp[9]  = cone[3]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4]; /* Front */
11851cfd6a4SMatthew G. Knepley         facesTmp[12] = cone[2]; facesTmp[13] = cone[1]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6]; /* Back */
11951cfd6a4SMatthew G. Knepley         facesTmp[16] = cone[3]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5]; /* Right */
12051cfd6a4SMatthew G. Knepley         facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[1]; /* Left */
1219a5b13a2SMatthew G Knepley         *faces = facesTmp;
1229a5b13a2SMatthew G Knepley       }
1239a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 6;
1249a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 4;
1259f074e33SMatthew G Knepley       break;
1269f074e33SMatthew G Knepley     default:
1279f074e33SMatthew G Knepley       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
1289f074e33SMatthew G Knepley     }
1299f074e33SMatthew G Knepley     break;
1309f074e33SMatthew G Knepley   default:
1319f074e33SMatthew G Knepley     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim);
1329f074e33SMatthew G Knepley   }
1339f074e33SMatthew G Knepley   PetscFunctionReturn(0);
1349f074e33SMatthew G Knepley }
1359f074e33SMatthew G Knepley 
1369a5b13a2SMatthew G Knepley /* This interpolates faces for cells at some stratum */
1379a5b13a2SMatthew G Knepley static PetscErrorCode DMPlexInterpolateFaces_Internal(DM dm, PetscInt cellDepth, DM idm)
1389f074e33SMatthew G Knepley {
139d4efc6f1SMatthew G. Knepley   DMLabel        subpointMap;
1409a5b13a2SMatthew G Knepley   PetscHashIJKL  faceTable;
1419a5b13a2SMatthew G Knepley   PetscInt      *pStart, *pEnd;
1429a5b13a2SMatthew G Knepley   PetscInt       cellDim, depth, faceDepth = cellDepth, numPoints = 0, faceSizeAll = 0, face, c, d;
1439f074e33SMatthew G Knepley   PetscErrorCode ierr;
1449f074e33SMatthew G Knepley 
1459f074e33SMatthew G Knepley   PetscFunctionBegin;
146c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &cellDim);CHKERRQ(ierr);
147d4efc6f1SMatthew G. Knepley   /* HACK: I need a better way to determine face dimension, or an alternative to GetFaces() */
148d4efc6f1SMatthew G. Knepley   ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr);
149d4efc6f1SMatthew G. Knepley   if (subpointMap) ++cellDim;
1509a5b13a2SMatthew G Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1519a5b13a2SMatthew G Knepley   ++depth;
1529a5b13a2SMatthew G Knepley   ++cellDepth;
1539a5b13a2SMatthew G Knepley   cellDim -= depth - cellDepth;
154dcca6d9dSJed Brown   ierr = PetscMalloc2(depth+1,&pStart,depth+1,&pEnd);CHKERRQ(ierr);
1559a5b13a2SMatthew G Knepley   for (d = depth-1; d >= faceDepth; --d) {
1569a5b13a2SMatthew G Knepley     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d+1], &pEnd[d+1]);CHKERRQ(ierr);
1579f074e33SMatthew G Knepley   }
1589a5b13a2SMatthew G Knepley   ierr = DMPlexGetDepthStratum(dm, -1, NULL, &pStart[faceDepth]);CHKERRQ(ierr);
1599a5b13a2SMatthew G Knepley   pEnd[faceDepth] = pStart[faceDepth];
1609a5b13a2SMatthew G Knepley   for (d = faceDepth-1; d >= 0; --d) {
1619a5b13a2SMatthew G Knepley     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
1629f074e33SMatthew G Knepley   }
1639a5b13a2SMatthew G Knepley   if (pEnd[cellDepth] > pStart[cellDepth]) {ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr);}
1649a5b13a2SMatthew 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);
1659a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr);
1669a5b13a2SMatthew G Knepley   for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) {
1679f074e33SMatthew G Knepley     const PetscInt *cellFaces;
168735a0255SMatthew G. Knepley     PetscInt        numCellFaces, faceSize, cf;
1699f074e33SMatthew G Knepley 
1709a5b13a2SMatthew G Knepley     ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
1719a5b13a2SMatthew G Knepley     if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll);
1729f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
1739a5b13a2SMatthew G Knepley       const PetscInt   *cellFace = &cellFaces[cf*faceSize];
1749a5b13a2SMatthew G Knepley       PetscHashIJKLKey  key;
175735a0255SMatthew G. Knepley       PetscHashIJKLIter missing, iter;
1769f074e33SMatthew G Knepley 
1779a5b13a2SMatthew G Knepley       if (faceSize == 2) {
1789a5b13a2SMatthew G Knepley         key.i = PetscMin(cellFace[0], cellFace[1]);
1799a5b13a2SMatthew G Knepley         key.j = PetscMax(cellFace[0], cellFace[1]);
1807fcd4ef0SJed Brown         key.k = 0;
1817fcd4ef0SJed Brown         key.l = 0;
1829a5b13a2SMatthew G Knepley       } else {
1839a5b13a2SMatthew G Knepley         key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0;
184302440fdSBarry Smith         ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr);
1859f074e33SMatthew G Knepley       }
186735a0255SMatthew G. Knepley       ierr = PetscHashIJKLPut(faceTable, key, &missing, &iter);CHKERRQ(ierr);
187735a0255SMatthew G. Knepley       if (missing) {ierr = PetscHashIJKLSet(faceTable, iter, face++);CHKERRQ(ierr);}
1889a5b13a2SMatthew G Knepley     }
189439ece16SMatthew G. Knepley     ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
1909a5b13a2SMatthew G Knepley   }
1919a5b13a2SMatthew G Knepley   pEnd[faceDepth] = face;
1929a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr);
1939a5b13a2SMatthew G Knepley   /* Count new points */
1949a5b13a2SMatthew G Knepley   for (d = 0; d <= depth; ++d) {
1959a5b13a2SMatthew G Knepley     numPoints += pEnd[d]-pStart[d];
1969a5b13a2SMatthew G Knepley   }
1979a5b13a2SMatthew G Knepley   ierr = DMPlexSetChart(idm, 0, numPoints);CHKERRQ(ierr);
1989a5b13a2SMatthew G Knepley   /* Set cone sizes */
1999a5b13a2SMatthew G Knepley   for (d = 0; d <= depth; ++d) {
2009a5b13a2SMatthew G Knepley     PetscInt coneSize, p;
2019f074e33SMatthew G Knepley 
2029a5b13a2SMatthew G Knepley     if (d == faceDepth) {
2039a5b13a2SMatthew G Knepley       for (p = pStart[d]; p < pEnd[d]; ++p) {
2049a5b13a2SMatthew G Knepley         /* I see no way to do this if we admit faces of different shapes */
2059a5b13a2SMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr);
2069a5b13a2SMatthew G Knepley       }
207a014044eSMatthew G Knepley     } else if (d == cellDepth) {
208a014044eSMatthew G Knepley       for (p = pStart[d]; p < pEnd[d]; ++p) {
209a014044eSMatthew G Knepley         /* Number of cell faces may be different from number of cell vertices*/
210a014044eSMatthew G Knepley         ierr = DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);CHKERRQ(ierr);
211a014044eSMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr);
212a014044eSMatthew G Knepley       }
2139a5b13a2SMatthew G Knepley     } else {
2149a5b13a2SMatthew G Knepley       for (p = pStart[d]; p < pEnd[d]; ++p) {
2159a5b13a2SMatthew G Knepley         ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
2169a5b13a2SMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr);
2179f074e33SMatthew G Knepley       }
2189f074e33SMatthew G Knepley     }
2199f074e33SMatthew G Knepley   }
2209f074e33SMatthew G Knepley   ierr = DMSetUp(idm);CHKERRQ(ierr);
2219f074e33SMatthew G Knepley   /* Get face cones from subsets of cell vertices */
2229a5b13a2SMatthew 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);
2239a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr);
2249a5b13a2SMatthew G Knepley   for (d = depth; d > cellDepth; --d) {
2259a5b13a2SMatthew G Knepley     const PetscInt *cone;
2269a5b13a2SMatthew G Knepley     PetscInt        p;
2279a5b13a2SMatthew G Knepley 
2289a5b13a2SMatthew G Knepley     for (p = pStart[d]; p < pEnd[d]; ++p) {
2299a5b13a2SMatthew G Knepley       ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
2309a5b13a2SMatthew G Knepley       ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr);
2319a5b13a2SMatthew G Knepley       ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr);
2329a5b13a2SMatthew G Knepley       ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr);
2339f074e33SMatthew G Knepley     }
2349a5b13a2SMatthew G Knepley   }
2359a5b13a2SMatthew G Knepley   for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) {
2369f074e33SMatthew G Knepley     const PetscInt *cellFaces;
237735a0255SMatthew G. Knepley     PetscInt        numCellFaces, faceSize, cf;
2389f074e33SMatthew G Knepley 
2399a5b13a2SMatthew G Knepley     ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
2409a5b13a2SMatthew G Knepley     if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll);
2419f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
2429a5b13a2SMatthew G Knepley       const PetscInt  *cellFace = &cellFaces[cf*faceSize];
2439a5b13a2SMatthew G Knepley       PetscHashIJKLKey key;
244735a0255SMatthew G. Knepley       PetscHashIJKLIter missing, iter;
2459f074e33SMatthew G Knepley 
2469a5b13a2SMatthew G Knepley       if (faceSize == 2) {
2479a5b13a2SMatthew G Knepley         key.i = PetscMin(cellFace[0], cellFace[1]);
2489a5b13a2SMatthew G Knepley         key.j = PetscMax(cellFace[0], cellFace[1]);
2497fcd4ef0SJed Brown         key.k = 0;
2507fcd4ef0SJed Brown         key.l = 0;
2519a5b13a2SMatthew G Knepley       } else {
2529a5b13a2SMatthew G Knepley         key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0;
253302440fdSBarry Smith         ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr);
2549f074e33SMatthew G Knepley       }
255735a0255SMatthew G. Knepley       ierr = PetscHashIJKLPut(faceTable, key, &missing, &iter);CHKERRQ(ierr);
256735a0255SMatthew G. Knepley       if (missing) {
2579a5b13a2SMatthew G Knepley         ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr);
258735a0255SMatthew G. Knepley         ierr = PetscHashIJKLSet(faceTable, iter, face);CHKERRQ(ierr);
259735a0255SMatthew G. Knepley         ierr = DMPlexInsertCone(idm, c, cf, face++);CHKERRQ(ierr);
2609a5b13a2SMatthew G Knepley       } else {
2619a5b13a2SMatthew G Knepley         const PetscInt *cone;
262735a0255SMatthew G. Knepley         PetscInt        coneSize, ornt, i, j, f;
2639f074e33SMatthew G Knepley 
264735a0255SMatthew G. Knepley         ierr = PetscHashIJKLGet(faceTable, iter, &f);CHKERRQ(ierr);
2659a5b13a2SMatthew G Knepley         ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr);
2662e1b13c2SMatthew G. Knepley         /* Orient face: Do not allow reverse orientation at the first vertex */
2679f074e33SMatthew G Knepley         ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr);
2689f074e33SMatthew G Knepley         ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr);
2699a5b13a2SMatthew 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);
2709a5b13a2SMatthew G Knepley         /* - First find the initial vertex */
2719a5b13a2SMatthew G Knepley         for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break;
2729a5b13a2SMatthew G Knepley         /* - Try forward comparison */
2739a5b13a2SMatthew G Knepley         for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break;
2749a5b13a2SMatthew G Knepley         if (j == faceSize) {
2759a5b13a2SMatthew G Knepley           if ((faceSize == 2) && (i == 1)) ornt = -2;
2769a5b13a2SMatthew G Knepley           else                             ornt = i;
2779a5b13a2SMatthew G Knepley         } else {
2789a5b13a2SMatthew G Knepley           /* - Try backward comparison */
2799a5b13a2SMatthew G Knepley           for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break;
2802e1b13c2SMatthew G. Knepley           if (j == faceSize) {
2812e1b13c2SMatthew G. Knepley             if (i == 0) ornt = -faceSize;
282dc1a705cSMatthew G. Knepley             else        ornt = -i;
2832e1b13c2SMatthew G. Knepley           } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation");
2849a5b13a2SMatthew G Knepley         }
2859a5b13a2SMatthew G Knepley         ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr);
2869f074e33SMatthew G Knepley       }
2879f074e33SMatthew G Knepley     }
288439ece16SMatthew G. Knepley     ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
2899f074e33SMatthew G Knepley   }
2909a5b13a2SMatthew 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]);
291c907b753SJed Brown   ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr);
2929a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr);
2936551a8c7SMatthew G. Knepley   ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr);
2949a5b13a2SMatthew G Knepley   ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr);
2959a5b13a2SMatthew G Knepley   ierr = DMPlexStratify(idm);CHKERRQ(ierr);
2969f074e33SMatthew G Knepley   PetscFunctionReturn(0);
2979f074e33SMatthew G Knepley }
2989f074e33SMatthew G Knepley 
2997bffde78SMichael Lange /* This interpolates the PointSF in parallel following local interpolation */
3007bffde78SMichael Lange static PetscErrorCode DMPlexInterpolatePointSF(DM dm, PetscSF pointSF, PetscInt depth)
3017bffde78SMichael Lange {
3029852e123SBarry Smith   PetscMPIInt        size, rank;
3037bffde78SMichael Lange   PetscInt           p, c, d, dof, offset;
304ffd6914dSSatish Balay   PetscInt           numLeaves, numRoots, candidatesSize, candidatesRemoteSize;
3057bffde78SMichael Lange   const PetscInt    *localPoints;
3067bffde78SMichael Lange   const PetscSFNode *remotePoints;
3077bffde78SMichael Lange   PetscSFNode       *candidates, *candidatesRemote, *claims;
3087bffde78SMichael Lange   PetscSection       candidateSection, candidateSectionRemote, claimSection;
3097bffde78SMichael Lange   PetscHashI         leafhash;
3107bffde78SMichael Lange   PetscHashIJ        roothash;
3117bffde78SMichael Lange   PetscHashIJKey     key;
3127bffde78SMichael Lange   PetscErrorCode     ierr;
3137bffde78SMichael Lange 
3147bffde78SMichael Lange   PetscFunctionBegin;
3159852e123SBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr);
3167bffde78SMichael Lange   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
3177bffde78SMichael Lange   ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr);
3189852e123SBarry Smith   if (size < 2 || numRoots < 0) PetscFunctionReturn(0);
31925afeb17SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr);
3207bffde78SMichael Lange   /* Build hashes of points in the SF for efficient lookup */
3217bffde78SMichael Lange   PetscHashICreate(leafhash);
3227bffde78SMichael Lange   PetscHashIJCreate(&roothash);
3237bffde78SMichael Lange   ierr = PetscHashIJSetMultivalued(roothash, PETSC_FALSE);CHKERRQ(ierr);
3247bffde78SMichael Lange   for (p = 0; p < numLeaves; ++p) {
3257bffde78SMichael Lange     PetscHashIAdd(leafhash, localPoints[p], p);
3267bffde78SMichael Lange     key.i = remotePoints[p].index; key.j = remotePoints[p].rank;
3277bffde78SMichael Lange     PetscHashIJAdd(roothash, key, p);
3287bffde78SMichael Lange   }
3297bffde78SMichael Lange   /* Build a section / SFNode array of candidate points in the single-level adjacency of leaves,
3307bffde78SMichael Lange      where each candidate is defined by the root entry for the other vertex that defines the edge. */
3317bffde78SMichael Lange   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSection);CHKERRQ(ierr);
3327bffde78SMichael Lange   ierr = PetscSectionSetChart(candidateSection, 0, numRoots);CHKERRQ(ierr);
3337bffde78SMichael Lange   {
334ffd6914dSSatish Balay     PetscInt leaf, root, idx, a, *adj = NULL;
3357bffde78SMichael Lange     for (p = 0; p < numLeaves; ++p) {
3367bffde78SMichael Lange       PetscInt adjSize = PETSC_DETERMINE;
3377bffde78SMichael Lange       ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr);
3387bffde78SMichael Lange       for (a = 0; a < adjSize; ++a) {
3397bffde78SMichael Lange         PetscHashIMap(leafhash, adj[a], leaf);
3407bffde78SMichael Lange         if (leaf >= 0) {ierr = PetscSectionAddDof(candidateSection, localPoints[p], 1);CHKERRQ(ierr);}
3417bffde78SMichael Lange       }
3427bffde78SMichael Lange     }
3437bffde78SMichael Lange     ierr = PetscSectionSetUp(candidateSection);CHKERRQ(ierr);
3447bffde78SMichael Lange     ierr = PetscSectionGetStorageSize(candidateSection, &candidatesSize);CHKERRQ(ierr);
3457bffde78SMichael Lange     ierr = PetscMalloc1(candidatesSize, &candidates);CHKERRQ(ierr);
3467bffde78SMichael Lange     for (p = 0; p < numLeaves; ++p) {
3477bffde78SMichael Lange       PetscInt adjSize = PETSC_DETERMINE;
3487bffde78SMichael Lange       ierr = PetscSectionGetOffset(candidateSection, localPoints[p], &offset);CHKERRQ(ierr);
3497bffde78SMichael Lange       ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr);
3507bffde78SMichael Lange       for (idx = 0, a = 0; a < adjSize; ++a) {
3517bffde78SMichael Lange         PetscHashIMap(leafhash, adj[a], root);
3527bffde78SMichael Lange         if (root >= 0) candidates[offset+idx++] = remotePoints[root];
3537bffde78SMichael Lange       }
3547bffde78SMichael Lange     }
3557bffde78SMichael Lange     ierr = PetscFree(adj);CHKERRQ(ierr);
3567bffde78SMichael Lange   }
3577bffde78SMichael Lange   /* Gather candidate section / array pair into the root partition via inverse(multi(pointSF)). */
3587bffde78SMichael Lange   {
3597bffde78SMichael Lange     PetscSF   sfMulti, sfInverse, sfCandidates;
3607bffde78SMichael Lange     PetscInt *remoteOffsets;
3617bffde78SMichael Lange     ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr);
3627bffde78SMichael Lange     ierr = PetscSFCreateInverseSF(sfMulti, &sfInverse);CHKERRQ(ierr);
3637bffde78SMichael Lange     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSectionRemote);CHKERRQ(ierr);
3647bffde78SMichael Lange     ierr = PetscSFDistributeSection(sfInverse, candidateSection, &remoteOffsets, candidateSectionRemote);CHKERRQ(ierr);
3657bffde78SMichael Lange     ierr = PetscSFCreateSectionSF(sfInverse, candidateSection, remoteOffsets, candidateSectionRemote, &sfCandidates);CHKERRQ(ierr);
3667bffde78SMichael Lange     ierr = PetscSectionGetStorageSize(candidateSectionRemote, &candidatesRemoteSize);CHKERRQ(ierr);
3677bffde78SMichael Lange     ierr = PetscMalloc1(candidatesRemoteSize, &candidatesRemote);CHKERRQ(ierr);
3687bffde78SMichael Lange     ierr = PetscSFBcastBegin(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr);
3697bffde78SMichael Lange     ierr = PetscSFBcastEnd(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr);
3707bffde78SMichael Lange     ierr = PetscSFDestroy(&sfInverse);CHKERRQ(ierr);
3717bffde78SMichael Lange     ierr = PetscSFDestroy(&sfCandidates);CHKERRQ(ierr);
3727bffde78SMichael Lange     ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
3737bffde78SMichael Lange   }
3747bffde78SMichael Lange   /* Walk local roots and check for each remote candidate whether we know all required points,
3757bffde78SMichael Lange      either from owning it or having a root entry in the point SF. If we do we place a claim
3767bffde78SMichael Lange      by replacing the vertex number with our edge ID. */
3777bffde78SMichael Lange   {
3787bffde78SMichael Lange     PetscInt        idx, root, joinSize, vertices[2];
3797bffde78SMichael Lange     const PetscInt *rootdegree, *join = NULL;
3807bffde78SMichael Lange     ierr = PetscSFComputeDegreeBegin(pointSF, &rootdegree);CHKERRQ(ierr);
3817bffde78SMichael Lange     ierr = PetscSFComputeDegreeEnd(pointSF, &rootdegree);CHKERRQ(ierr);
3827bffde78SMichael Lange     /* Loop remote edge connections and put in a claim if both vertices are known */
3837bffde78SMichael Lange     for (idx = 0, p = 0; p < numRoots; ++p) {
3847bffde78SMichael Lange       for (d = 0; d < rootdegree[p]; ++d) {
3857bffde78SMichael Lange         ierr = PetscSectionGetDof(candidateSectionRemote, idx, &dof);CHKERRQ(ierr);
3867bffde78SMichael Lange         ierr = PetscSectionGetOffset(candidateSectionRemote, idx, &offset);CHKERRQ(ierr);
3877bffde78SMichael Lange         for (c = 0; c < dof; ++c) {
3887bffde78SMichael Lange           /* We own both vertices, so we claim the edge by replacing vertex with edge */
3897bffde78SMichael Lange           if (candidatesRemote[offset+c].rank == rank) {
3907bffde78SMichael Lange             vertices[0] = p; vertices[1] = candidatesRemote[offset+c].index;
3917bffde78SMichael Lange             ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
3927bffde78SMichael Lange             if (joinSize == 1) candidatesRemote[offset+c].index = join[0];
3937bffde78SMichael Lange             ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
3947bffde78SMichael Lange             continue;
3957bffde78SMichael Lange           }
3967bffde78SMichael Lange           /* If we own one vertex and share a root with the other, we claim it */
3977bffde78SMichael Lange           key.i = candidatesRemote[offset+c].index; key.j = candidatesRemote[offset+c].rank;
3987bffde78SMichael Lange           PetscHashIJGet(roothash, key, &root);
3997bffde78SMichael Lange           if (root >= 0) {
4007bffde78SMichael Lange             vertices[0] = p; vertices[1] = localPoints[root];
4017bffde78SMichael Lange             ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
4027bffde78SMichael Lange             if (joinSize == 1) {
4037bffde78SMichael Lange               candidatesRemote[offset+c].index = join[0];
4047bffde78SMichael Lange               candidatesRemote[offset+c].rank = rank;
4057bffde78SMichael Lange             }
4067bffde78SMichael Lange             ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
4077bffde78SMichael Lange           }
4087bffde78SMichael Lange         }
4097bffde78SMichael Lange         idx++;
4107bffde78SMichael Lange       }
4117bffde78SMichael Lange     }
4127bffde78SMichael Lange   }
4137bffde78SMichael Lange   /* Push claims back to receiver via the MultiSF and derive new pointSF mapping on receiver */
4147bffde78SMichael Lange   {
4157bffde78SMichael Lange     PetscSF         sfMulti, sfClaims, sfPointNew;
4167bffde78SMichael Lange     PetscHashI      claimshash;
4177bffde78SMichael Lange     PetscInt        size, pStart, pEnd, root, joinSize, numLocalNew;
4187bffde78SMichael Lange     PetscInt       *remoteOffsets, *localPointsNew, vertices[2];
419ffd6914dSSatish Balay     const PetscInt *join = NULL;
4207bffde78SMichael Lange     PetscSFNode    *remotePointsNew;
4217bffde78SMichael Lange     ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr);
4227bffde78SMichael Lange     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &claimSection);CHKERRQ(ierr);
4237bffde78SMichael Lange     ierr = PetscSFDistributeSection(sfMulti, candidateSectionRemote, &remoteOffsets, claimSection);CHKERRQ(ierr);
4247bffde78SMichael Lange     ierr = PetscSFCreateSectionSF(sfMulti, candidateSectionRemote, remoteOffsets, claimSection, &sfClaims);CHKERRQ(ierr);
4257bffde78SMichael Lange     ierr = PetscSectionGetStorageSize(claimSection, &size);CHKERRQ(ierr);
4267bffde78SMichael Lange     ierr = PetscMalloc1(size, &claims);CHKERRQ(ierr);
4277bffde78SMichael Lange     ierr = PetscSFBcastBegin(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr);
4287bffde78SMichael Lange     ierr = PetscSFBcastEnd(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr);
4297bffde78SMichael Lange     ierr = PetscSFDestroy(&sfClaims);CHKERRQ(ierr);
4307bffde78SMichael Lange     ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
4317bffde78SMichael Lange     /* Walk the original section of local supports and add an SF entry for each updated item */
4327bffde78SMichael Lange     PetscHashICreate(claimshash);
4337bffde78SMichael Lange     for (p = 0; p < numRoots; ++p) {
4347bffde78SMichael Lange       ierr = PetscSectionGetDof(candidateSection, p, &dof);CHKERRQ(ierr);
4357bffde78SMichael Lange       ierr = PetscSectionGetOffset(candidateSection, p, &offset);CHKERRQ(ierr);
4367bffde78SMichael Lange       for (d = 0; d < dof; ++d) {
4377bffde78SMichael Lange         if (candidates[offset+d].index != claims[offset+d].index) {
4387bffde78SMichael Lange           key.i = candidates[offset+d].index; key.j = candidates[offset+d].rank;
4397bffde78SMichael Lange           PetscHashIJGet(roothash, key, &root);
4407bffde78SMichael Lange           if (root >= 0) {
4417bffde78SMichael Lange             vertices[0] = p; vertices[1] = localPoints[root];
4427bffde78SMichael Lange             ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
4437bffde78SMichael Lange             if (joinSize == 1) PetscHashIAdd(claimshash, join[0], offset+d);
4447bffde78SMichael Lange             ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
4457bffde78SMichael Lange           }
4467bffde78SMichael Lange         }
4477bffde78SMichael Lange       }
4487bffde78SMichael Lange     }
4497bffde78SMichael Lange     /* Create new pointSF from hashed claims */
4507bffde78SMichael Lange     PetscHashISize(claimshash, numLocalNew);
4517bffde78SMichael Lange     ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
4527bffde78SMichael Lange     ierr = PetscMalloc1(numLeaves + numLocalNew, &localPointsNew);CHKERRQ(ierr);
4537bffde78SMichael Lange     ierr = PetscMalloc1(numLeaves + numLocalNew, &remotePointsNew);CHKERRQ(ierr);
4547bffde78SMichael Lange     for (p = 0; p < numLeaves; ++p) {
4557bffde78SMichael Lange       localPointsNew[p] = localPoints[p];
4567bffde78SMichael Lange       remotePointsNew[p].index = remotePoints[p].index;
4577bffde78SMichael Lange       remotePointsNew[p].rank = remotePoints[p].rank;
4587bffde78SMichael Lange     }
459f3190fc2SToby Isaac     p = numLeaves;
460f3190fc2SToby Isaac     ierr = PetscHashIGetKeys(claimshash, &p, localPointsNew);CHKERRQ(ierr);
461f3190fc2SToby Isaac     ierr = PetscSortInt(numLocalNew,&localPointsNew[numLeaves]);CHKERRQ(ierr);
4627bffde78SMichael Lange     for (p = numLeaves; p < numLeaves + numLocalNew; ++p) {
4637bffde78SMichael Lange       PetscHashIMap(claimshash, localPointsNew[p], offset);
4647bffde78SMichael Lange       remotePointsNew[p] = claims[offset];
4657bffde78SMichael Lange     }
4667bffde78SMichael Lange     ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), &sfPointNew);CHKERRQ(ierr);
4677bffde78SMichael Lange     ierr = PetscSFSetGraph(sfPointNew, pEnd-pStart, numLeaves+numLocalNew, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr);
4687bffde78SMichael Lange     ierr = DMSetPointSF(dm, sfPointNew);CHKERRQ(ierr);
4697bffde78SMichael Lange     ierr = PetscSFDestroy(&sfPointNew);CHKERRQ(ierr);
4707bffde78SMichael Lange     PetscHashIDestroy(claimshash);
4717bffde78SMichael Lange   }
4727bffde78SMichael Lange   PetscHashIDestroy(leafhash);
4737bffde78SMichael Lange   ierr = PetscHashIJDestroy(&roothash);CHKERRQ(ierr);
4747bffde78SMichael Lange   ierr = PetscSectionDestroy(&candidateSection);CHKERRQ(ierr);
4757bffde78SMichael Lange   ierr = PetscSectionDestroy(&candidateSectionRemote);CHKERRQ(ierr);
4767bffde78SMichael Lange   ierr = PetscSectionDestroy(&claimSection);CHKERRQ(ierr);
4777bffde78SMichael Lange   ierr = PetscFree(candidates);CHKERRQ(ierr);
4787bffde78SMichael Lange   ierr = PetscFree(candidatesRemote);CHKERRQ(ierr);
4797bffde78SMichael Lange   ierr = PetscFree(claims);CHKERRQ(ierr);
48025afeb17SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr);
4817bffde78SMichael Lange   PetscFunctionReturn(0);
4827bffde78SMichael Lange }
4837bffde78SMichael Lange 
4840c795ddcSMatthew G. Knepley /*@C
48580330477SMatthew G. Knepley   DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc.
48680330477SMatthew G. Knepley 
48780330477SMatthew G. Knepley   Collective on DM
48880330477SMatthew G. Knepley 
489e56d480eSMatthew G. Knepley   Input Parameters:
490e56d480eSMatthew G. Knepley + dm - The DMPlex object with only cells and vertices
49110a67516SMatthew G. Knepley - dmInt - The interpolated DM
49280330477SMatthew G. Knepley 
49380330477SMatthew G. Knepley   Output Parameter:
4944e3744c5SMatthew G. Knepley . dmInt - The complete DMPlex object
49580330477SMatthew G. Knepley 
49680330477SMatthew G. Knepley   Level: intermediate
49780330477SMatthew G. Knepley 
49880330477SMatthew G. Knepley .keywords: mesh
4994e3744c5SMatthew G. Knepley .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList()
50080330477SMatthew G. Knepley @*/
5019f074e33SMatthew G Knepley PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt)
5029f074e33SMatthew G Knepley {
5039a5b13a2SMatthew G Knepley   DM             idm, odm = dm;
5047bffde78SMichael Lange   PetscSF        sfPoint;
5052e1b13c2SMatthew G. Knepley   PetscInt       depth, dim, d;
50610a67516SMatthew G. Knepley   const char    *name;
5079f074e33SMatthew G Knepley   PetscErrorCode ierr;
5089f074e33SMatthew G Knepley 
5099f074e33SMatthew G Knepley   PetscFunctionBegin;
51010a67516SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
51110a67516SMatthew G. Knepley   PetscValidPointer(dmInt, 2);
512a72f3261SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr);
5132e1b13c2SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
514c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
51576b791aaSMatthew G. Knepley   if (dim <= 1) {
51676b791aaSMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
51776b791aaSMatthew G. Knepley     idm  = dm;
51876b791aaSMatthew G. Knepley   }
5199a5b13a2SMatthew G Knepley   for (d = 1; d < dim; ++d) {
5209a5b13a2SMatthew G Knepley     /* Create interpolated mesh */
52110a67516SMatthew G. Knepley     ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr);
5229a5b13a2SMatthew G Knepley     ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr);
523c73cfb54SMatthew G. Knepley     ierr = DMSetDimension(idm, dim);CHKERRQ(ierr);
5247bffde78SMichael Lange     if (depth > 0) {
5257bffde78SMichael Lange       ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr);
5267bffde78SMichael Lange       ierr = DMGetPointSF(odm, &sfPoint);CHKERRQ(ierr);
5277bffde78SMichael Lange       ierr = DMPlexInterpolatePointSF(idm, sfPoint, depth);CHKERRQ(ierr);
5287bffde78SMichael Lange     }
5299a5b13a2SMatthew G Knepley     if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);}
5309a5b13a2SMatthew G Knepley     odm  = idm;
5319f074e33SMatthew G Knepley   }
53210a67516SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) dm,  &name);CHKERRQ(ierr);
53310a67516SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) idm,  name);CHKERRQ(ierr);
53410a67516SMatthew G. Knepley   ierr = DMPlexCopyCoordinates(dm, idm);CHKERRQ(ierr);
53510a67516SMatthew G. Knepley   ierr = DMCopyLabels(dm, idm);CHKERRQ(ierr);
5369a5b13a2SMatthew G Knepley   *dmInt = idm;
537a72f3261SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr);
5389f074e33SMatthew G Knepley   PetscFunctionReturn(0);
5399f074e33SMatthew G Knepley }
54007b0a1fcSMatthew G Knepley 
54180330477SMatthew G. Knepley /*@
54280330477SMatthew G. Knepley   DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices
54380330477SMatthew G. Knepley 
54480330477SMatthew G. Knepley   Collective on DM
54580330477SMatthew G. Knepley 
54680330477SMatthew G. Knepley   Input Parameter:
54780330477SMatthew G. Knepley . dmA - The DMPlex object with initial coordinates
54880330477SMatthew G. Knepley 
54980330477SMatthew G. Knepley   Output Parameter:
55080330477SMatthew G. Knepley . dmB - The DMPlex object with copied coordinates
55180330477SMatthew G. Knepley 
55280330477SMatthew G. Knepley   Level: intermediate
55380330477SMatthew G. Knepley 
55480330477SMatthew G. Knepley   Note: This is typically used when adding pieces other than vertices to a mesh
55580330477SMatthew G. Knepley 
55680330477SMatthew G. Knepley .keywords: mesh
55765f90189SJed Brown .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMGetCoordinateSection()
55880330477SMatthew G. Knepley @*/
55907b0a1fcSMatthew G Knepley PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB)
56007b0a1fcSMatthew G Knepley {
56107b0a1fcSMatthew G Knepley   Vec            coordinatesA, coordinatesB;
56207b0a1fcSMatthew G Knepley   PetscSection   coordSectionA, coordSectionB;
56307b0a1fcSMatthew G Knepley   PetscScalar   *coordsA, *coordsB;
564*0bedd6beSMatthew G. Knepley   PetscInt       spaceDim, Nf, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d;
56507b0a1fcSMatthew G Knepley   PetscErrorCode ierr;
56607b0a1fcSMatthew G Knepley 
56707b0a1fcSMatthew G Knepley   PetscFunctionBegin;
56876b791aaSMatthew G. Knepley   if (dmA == dmB) PetscFunctionReturn(0);
56907b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr);
57007b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr);
57107b0a1fcSMatthew 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);
57269d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr);
57369d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr);
574972bc18aSToby Isaac   if (coordSectionA == coordSectionB) PetscFunctionReturn(0);
575*0bedd6beSMatthew G. Knepley   ierr = PetscSectionGetNumFields(coordSectionA, &Nf);CHKERRQ(ierr);
576*0bedd6beSMatthew G. Knepley   if (!Nf) PetscFunctionReturn(0);
577*0bedd6beSMatthew G. Knepley   if (Nf > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of coordinate fields must be 1, not %D", Nf);
578df26b574SMatthew G. Knepley   if (!coordSectionB) {
579df26b574SMatthew G. Knepley     PetscInt dim;
580df26b574SMatthew G. Knepley 
581df26b574SMatthew G. Knepley     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) coordSectionA), &coordSectionB);CHKERRQ(ierr);
582df26b574SMatthew G. Knepley     ierr = DMGetCoordinateDim(dmA, &dim);CHKERRQ(ierr);
583df26b574SMatthew G. Knepley     ierr = DMSetCoordinateSection(dmB, dim, coordSectionB);CHKERRQ(ierr);
584df26b574SMatthew G. Knepley     ierr = PetscObjectDereference((PetscObject) coordSectionB);CHKERRQ(ierr);
585df26b574SMatthew G. Knepley   }
58607b0a1fcSMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr);
58707b0a1fcSMatthew G Knepley   ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr);
58807b0a1fcSMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr);
58907b0a1fcSMatthew G Knepley   ierr = PetscSectionSetChart(coordSectionB, vStartB, vEndB);CHKERRQ(ierr);
59007b0a1fcSMatthew G Knepley   for (v = vStartB; v < vEndB; ++v) {
59107b0a1fcSMatthew G Knepley     ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr);
59207b0a1fcSMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr);
59307b0a1fcSMatthew G Knepley   }
59407b0a1fcSMatthew G Knepley   ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr);
59507b0a1fcSMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr);
59607b0a1fcSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr);
5978b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinatesB);CHKERRQ(ierr);
59807b0a1fcSMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr);
59907b0a1fcSMatthew G Knepley   ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr);
6002eb5907fSJed Brown   ierr = VecSetType(coordinatesB,VECSTANDARD);CHKERRQ(ierr);
60107b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr);
60207b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr);
60307b0a1fcSMatthew G Knepley   for (v = 0; v < vEndB-vStartB; ++v) {
60407b0a1fcSMatthew G Knepley     for (d = 0; d < spaceDim; ++d) {
60507b0a1fcSMatthew G Knepley       coordsB[v*spaceDim+d] = coordsA[v*spaceDim+d];
60607b0a1fcSMatthew G Knepley     }
60707b0a1fcSMatthew G Knepley   }
60807b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr);
60907b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr);
61007b0a1fcSMatthew G Knepley   ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr);
61107b0a1fcSMatthew G Knepley   ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr);
61207b0a1fcSMatthew G Knepley   PetscFunctionReturn(0);
61307b0a1fcSMatthew G Knepley }
6145c386225SMatthew G. Knepley 
6154e3744c5SMatthew G. Knepley /*@
6164e3744c5SMatthew G. Knepley   DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh
6174e3744c5SMatthew G. Knepley 
6184e3744c5SMatthew G. Knepley   Collective on DM
6194e3744c5SMatthew G. Knepley 
6204e3744c5SMatthew G. Knepley   Input Parameter:
6214e3744c5SMatthew G. Knepley . dm - The complete DMPlex object
6224e3744c5SMatthew G. Knepley 
6234e3744c5SMatthew G. Knepley   Output Parameter:
6244e3744c5SMatthew G. Knepley . dmUnint - The DMPlex object with only cells and vertices
6254e3744c5SMatthew G. Knepley 
6264e3744c5SMatthew G. Knepley   Level: intermediate
6274e3744c5SMatthew G. Knepley 
6284e3744c5SMatthew G. Knepley .keywords: mesh
6294e3744c5SMatthew G. Knepley .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList()
6304e3744c5SMatthew G. Knepley @*/
6314e3744c5SMatthew G. Knepley PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint)
6324e3744c5SMatthew G. Knepley {
6334e3744c5SMatthew G. Knepley   DM             udm;
634595d4abbSMatthew G. Knepley   PetscInt       dim, vStart, vEnd, cStart, cEnd, c, maxConeSize = 0, *cone;
6354e3744c5SMatthew G. Knepley   PetscErrorCode ierr;
6364e3744c5SMatthew G. Knepley 
6374e3744c5SMatthew G. Knepley   PetscFunctionBegin;
638c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
6394e3744c5SMatthew G. Knepley   if (dim <= 1) {
6404e3744c5SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
641595d4abbSMatthew G. Knepley     *dmUnint = dm;
642595d4abbSMatthew G. Knepley     PetscFunctionReturn(0);
6434e3744c5SMatthew G. Knepley   }
644595d4abbSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
645595d4abbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
6464e3744c5SMatthew G. Knepley   ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr);
6474e3744c5SMatthew G. Knepley   ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr);
648c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(udm, dim);CHKERRQ(ierr);
6494e3744c5SMatthew G. Knepley   ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr);
6504e3744c5SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
651595d4abbSMatthew G. Knepley     PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
6524e3744c5SMatthew G. Knepley 
6534e3744c5SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
6544e3744c5SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
6554e3744c5SMatthew G. Knepley       const PetscInt p = closure[cl];
6564e3744c5SMatthew G. Knepley 
6574e3744c5SMatthew G. Knepley       if ((p >= vStart) && (p < vEnd)) ++coneSize;
6584e3744c5SMatthew G. Knepley     }
6594e3744c5SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
6604e3744c5SMatthew G. Knepley     ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr);
661595d4abbSMatthew G. Knepley     maxConeSize = PetscMax(maxConeSize, coneSize);
6624e3744c5SMatthew G. Knepley   }
6634e3744c5SMatthew G. Knepley   ierr = DMSetUp(udm);CHKERRQ(ierr);
664785e854fSJed Brown   ierr = PetscMalloc1(maxConeSize, &cone);CHKERRQ(ierr);
6654e3744c5SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
666595d4abbSMatthew G. Knepley     PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
6674e3744c5SMatthew G. Knepley 
6684e3744c5SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
6694e3744c5SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
6704e3744c5SMatthew G. Knepley       const PetscInt p = closure[cl];
6714e3744c5SMatthew G. Knepley 
6724e3744c5SMatthew G. Knepley       if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p;
6734e3744c5SMatthew G. Knepley     }
6744e3744c5SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
6754e3744c5SMatthew G. Knepley     ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr);
6764e3744c5SMatthew G. Knepley   }
6774e3744c5SMatthew G. Knepley   ierr = PetscFree(cone);CHKERRQ(ierr);
6784e3744c5SMatthew G. Knepley   ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr);
6794e3744c5SMatthew G. Knepley   ierr = DMPlexStratify(udm);CHKERRQ(ierr);
6805c7de58aSMatthew G. Knepley   /* Reduce SF */
6815c7de58aSMatthew G. Knepley   {
6825c7de58aSMatthew G. Knepley     PetscSF            sfPoint, sfPointUn;
6835c7de58aSMatthew G. Knepley     const PetscSFNode *remotePoints;
6845c7de58aSMatthew G. Knepley     const PetscInt    *localPoints;
6855c7de58aSMatthew G. Knepley     PetscSFNode       *remotePointsUn;
6865c7de58aSMatthew G. Knepley     PetscInt          *localPointsUn;
6875c7de58aSMatthew G. Knepley     PetscInt           vEnd, numRoots, numLeaves, l;
6885c7de58aSMatthew G. Knepley     PetscInt           numLeavesUn = 0, n = 0;
6895c7de58aSMatthew G. Knepley     PetscErrorCode     ierr;
6905c7de58aSMatthew G. Knepley 
6915c7de58aSMatthew G. Knepley     /* Get original SF information */
6925c7de58aSMatthew G. Knepley     ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
6935c7de58aSMatthew G. Knepley     ierr = DMGetPointSF(udm, &sfPointUn);CHKERRQ(ierr);
6945c7de58aSMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, 0, NULL, &vEnd);CHKERRQ(ierr);
6955c7de58aSMatthew G. Knepley     ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr);
6965c7de58aSMatthew G. Knepley     /* Allocate space for cells and vertices */
6975c7de58aSMatthew G. Knepley     for (l = 0; l < numLeaves; ++l) if (localPoints[l] < vEnd) numLeavesUn++;
6985c7de58aSMatthew G. Knepley     /* Fill in leaves */
6995c7de58aSMatthew G. Knepley     if (vEnd >= 0) {
7005c7de58aSMatthew G. Knepley       ierr = PetscMalloc1(numLeavesUn, &remotePointsUn);CHKERRQ(ierr);
7015c7de58aSMatthew G. Knepley       ierr = PetscMalloc1(numLeavesUn, &localPointsUn);CHKERRQ(ierr);
7025c7de58aSMatthew G. Knepley       for (l = 0; l < numLeaves; l++) {
7035c7de58aSMatthew G. Knepley         if (localPoints[l] < vEnd) {
7045c7de58aSMatthew G. Knepley           localPointsUn[n]        = localPoints[l];
7055c7de58aSMatthew G. Knepley           remotePointsUn[n].rank  = remotePoints[l].rank;
7065c7de58aSMatthew G. Knepley           remotePointsUn[n].index = remotePoints[l].index;
7075c7de58aSMatthew G. Knepley           ++n;
7085c7de58aSMatthew G. Knepley         }
7095c7de58aSMatthew G. Knepley       }
7105c7de58aSMatthew G. Knepley       if (n != numLeavesUn) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent number of leaves %d != %d", n, numLeavesUn);
7115c7de58aSMatthew G. Knepley       ierr = PetscSFSetGraph(sfPointUn, vEnd, numLeavesUn, localPointsUn, PETSC_OWN_POINTER, remotePointsUn, PETSC_OWN_POINTER);CHKERRQ(ierr);
7125c7de58aSMatthew G. Knepley     }
7135c7de58aSMatthew G. Knepley   }
7144e3744c5SMatthew G. Knepley   *dmUnint = udm;
7154e3744c5SMatthew G. Knepley   PetscFunctionReturn(0);
7164e3744c5SMatthew G. Knepley }
717