xref: /petsc/src/dm/impls/plex/plexinterpolate.c (revision 4e3744c5f6311de88db7a2dd30c820908180377f)
19f074e33SMatthew G Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
29f074e33SMatthew G Knepley #include <../src/sys/utils/hash.h>
39f074e33SMatthew G Knepley 
49f074e33SMatthew G Knepley #undef __FUNCT__
59a5b13a2SMatthew G Knepley #define __FUNCT__ "DMPlexGetFaces_Internal"
69f074e33SMatthew G Knepley /*
79a5b13a2SMatthew G Knepley   DMPlexGetFaces_Internal - Gets groups of vertices that correspond to faces for the given cell
89f074e33SMatthew G Knepley */
99a5b13a2SMatthew G Knepley static 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       *facesTmp;
139a5b13a2SMatthew G Knepley   PetscInt        maxConeSize, maxSupportSize, coneSize;
149f074e33SMatthew G Knepley   PetscErrorCode  ierr;
159f074e33SMatthew G Knepley 
169f074e33SMatthew G Knepley   PetscFunctionBegin;
179f074e33SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
189a5b13a2SMatthew G Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
199a5b13a2SMatthew G Knepley   ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), PETSC_INT, &facesTmp);CHKERRQ(ierr);
209f074e33SMatthew G Knepley   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
219f074e33SMatthew G Knepley   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
229f074e33SMatthew G Knepley   switch (dim) {
239f074e33SMatthew G Knepley   case 2:
249f074e33SMatthew G Knepley     switch (coneSize) {
259f074e33SMatthew G Knepley     case 3:
269a5b13a2SMatthew G Knepley       if (faces) {
279a5b13a2SMatthew G Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
289a5b13a2SMatthew G Knepley         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
299a5b13a2SMatthew G Knepley         facesTmp[4] = cone[2]; facesTmp[5] = cone[0];
309a5b13a2SMatthew G Knepley         *faces = facesTmp;
319a5b13a2SMatthew G Knepley       }
329a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 3;
339a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 2;
349f074e33SMatthew G Knepley       break;
359f074e33SMatthew G Knepley     case 4:
369a5b13a2SMatthew G Knepley       /* Vertices follow right hand rule */
379a5b13a2SMatthew G Knepley       if (faces) {
389a5b13a2SMatthew G Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
399a5b13a2SMatthew G Knepley         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
409a5b13a2SMatthew G Knepley         facesTmp[4] = cone[2]; facesTmp[5] = cone[3];
419a5b13a2SMatthew G Knepley         facesTmp[6] = cone[3]; facesTmp[7] = cone[0];
429a5b13a2SMatthew G Knepley         *faces = facesTmp;
439a5b13a2SMatthew G Knepley       }
449a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 4;
459a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 2;
469a5b13a2SMatthew G Knepley       if (faces)    *faces            = facesTmp;
479f074e33SMatthew G Knepley       break;
489f074e33SMatthew G Knepley     default:
499f074e33SMatthew G Knepley       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
509f074e33SMatthew G Knepley     }
519f074e33SMatthew G Knepley     break;
529f074e33SMatthew G Knepley   case 3:
539f074e33SMatthew G Knepley     switch (coneSize) {
549f074e33SMatthew G Knepley     case 3:
559a5b13a2SMatthew G Knepley       if (faces) {
569a5b13a2SMatthew G Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
579a5b13a2SMatthew G Knepley         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
589a5b13a2SMatthew G Knepley         facesTmp[4] = cone[2]; facesTmp[5] = cone[0];
599a5b13a2SMatthew G Knepley         *faces = facesTmp;
609a5b13a2SMatthew G Knepley       }
619a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 3;
629a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 2;
639a5b13a2SMatthew G Knepley       if (faces)    *faces            = facesTmp;
649f074e33SMatthew G Knepley       break;
659f074e33SMatthew G Knepley     case 4:
662e1b13c2SMatthew G. Knepley       /* Vertices of first face follow right hand rule and normal points away from last vertex */
679a5b13a2SMatthew G Knepley       if (faces) {
682e1b13c2SMatthew G. Knepley         facesTmp[0] = cone[0]; facesTmp[1]  = cone[1]; facesTmp[2]  = cone[2];
692e1b13c2SMatthew G. Knepley         facesTmp[3] = cone[0]; facesTmp[4]  = cone[3]; facesTmp[5]  = cone[1];
702e1b13c2SMatthew G. Knepley         facesTmp[6] = cone[0]; facesTmp[7]  = cone[2]; facesTmp[8]  = cone[3];
712e1b13c2SMatthew G. Knepley         facesTmp[9] = cone[2]; facesTmp[10] = cone[1]; facesTmp[11] = cone[3];
729a5b13a2SMatthew G Knepley         *faces = facesTmp;
739a5b13a2SMatthew G Knepley       }
749a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 4;
759a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 3;
769a5b13a2SMatthew G Knepley       if (faces)    *faces            = facesTmp;
779a5b13a2SMatthew G Knepley       break;
789a5b13a2SMatthew G Knepley     case 8:
799a5b13a2SMatthew G Knepley       if (faces) {
802e1b13c2SMatthew G. Knepley         facesTmp[0]  = cone[0]; facesTmp[1]  = cone[1]; facesTmp[2]  = cone[2]; facesTmp[3]  = cone[3];
81a014044eSMatthew G Knepley         facesTmp[4]  = cone[4]; facesTmp[5]  = cone[5]; facesTmp[6]  = cone[6]; facesTmp[7]  = cone[7];
822e1b13c2SMatthew G. Knepley         facesTmp[8]  = cone[0]; facesTmp[9]  = cone[3]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4];
832e1b13c2SMatthew G. Knepley         facesTmp[12] = cone[2]; facesTmp[13] = cone[1]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6];
842e1b13c2SMatthew G. Knepley         facesTmp[16] = cone[3]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5];
852e1b13c2SMatthew G. Knepley         facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[1];
869a5b13a2SMatthew G Knepley         *faces = facesTmp;
879a5b13a2SMatthew G Knepley       }
889a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 6;
899a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 4;
909f074e33SMatthew G Knepley       break;
919f074e33SMatthew G Knepley     default:
929f074e33SMatthew G Knepley       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
939f074e33SMatthew G Knepley     }
949f074e33SMatthew G Knepley     break;
959f074e33SMatthew G Knepley   default:
969f074e33SMatthew G Knepley     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim);
979f074e33SMatthew G Knepley   }
989a5b13a2SMatthew G Knepley   ierr = DMRestoreWorkArray(dm, 0, PETSC_INT, &facesTmp);CHKERRQ(ierr);
999f074e33SMatthew G Knepley   PetscFunctionReturn(0);
1009f074e33SMatthew G Knepley }
1019f074e33SMatthew G Knepley 
1029f074e33SMatthew G Knepley #undef __FUNCT__
1039a5b13a2SMatthew G Knepley #define __FUNCT__ "DMPlexInterpolateFaces_Internal"
1049a5b13a2SMatthew G Knepley /* This interpolates faces for cells at some stratum */
1059a5b13a2SMatthew G Knepley static PetscErrorCode DMPlexInterpolateFaces_Internal(DM dm, PetscInt cellDepth, DM idm)
1069f074e33SMatthew G Knepley {
107d4efc6f1SMatthew G. Knepley   DMLabel        subpointMap;
1089a5b13a2SMatthew G Knepley   PetscHashIJKL  faceTable;
1099a5b13a2SMatthew G Knepley   PetscInt      *pStart, *pEnd;
1109a5b13a2SMatthew G Knepley   PetscInt       cellDim, depth, faceDepth = cellDepth, numPoints = 0, faceSizeAll = 0, face, c, d;
1119f074e33SMatthew G Knepley   PetscErrorCode ierr;
1129f074e33SMatthew G Knepley 
1139f074e33SMatthew G Knepley   PetscFunctionBegin;
1149a5b13a2SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &cellDim);CHKERRQ(ierr);
115d4efc6f1SMatthew G. Knepley   /* HACK: I need a better way to determine face dimension, or an alternative to GetFaces() */
116d4efc6f1SMatthew G. Knepley   ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr);
117d4efc6f1SMatthew G. Knepley   if (subpointMap) ++cellDim;
1189a5b13a2SMatthew G Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1199a5b13a2SMatthew G Knepley   ++depth;
1209a5b13a2SMatthew G Knepley   ++cellDepth;
1219a5b13a2SMatthew G Knepley   cellDim -= depth - cellDepth;
1229a5b13a2SMatthew G Knepley   ierr = PetscMalloc2(depth+1,PetscInt,&pStart,depth+1,PetscInt,&pEnd);CHKERRQ(ierr);
1239a5b13a2SMatthew G Knepley   for (d = depth-1; d >= faceDepth; --d) {
1249a5b13a2SMatthew G Knepley     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d+1], &pEnd[d+1]);CHKERRQ(ierr);
1259f074e33SMatthew G Knepley   }
1269a5b13a2SMatthew G Knepley   ierr = DMPlexGetDepthStratum(dm, -1, NULL, &pStart[faceDepth]);CHKERRQ(ierr);
1279a5b13a2SMatthew G Knepley   pEnd[faceDepth] = pStart[faceDepth];
1289a5b13a2SMatthew G Knepley   for (d = faceDepth-1; d >= 0; --d) {
1299a5b13a2SMatthew G Knepley     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
1309f074e33SMatthew G Knepley   }
1319a5b13a2SMatthew G Knepley   if (pEnd[cellDepth] > pStart[cellDepth]) {ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr);}
1329a5b13a2SMatthew 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);
1339a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr);
1349a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLSetMultivalued(faceTable, PETSC_FALSE);CHKERRQ(ierr);
1359a5b13a2SMatthew G Knepley   for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) {
1369f074e33SMatthew G Knepley     const PetscInt *cellFaces;
1379a5b13a2SMatthew G Knepley     PetscInt        numCellFaces, faceSize, cf, f;
1389f074e33SMatthew G Knepley 
1399a5b13a2SMatthew G Knepley     ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
1409a5b13a2SMatthew G Knepley     if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll);
1419f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
1429a5b13a2SMatthew G Knepley       const PetscInt  *cellFace = &cellFaces[cf*faceSize];
1439a5b13a2SMatthew G Knepley       PetscHashIJKLKey key;
1449f074e33SMatthew G Knepley 
1459a5b13a2SMatthew G Knepley       if (faceSize == 2) {
1469a5b13a2SMatthew G Knepley         key.i = PetscMin(cellFace[0], cellFace[1]);
1479a5b13a2SMatthew G Knepley         key.j = PetscMax(cellFace[0], cellFace[1]);
1487fcd4ef0SJed Brown         key.k = 0;
1497fcd4ef0SJed Brown         key.l = 0;
1509a5b13a2SMatthew G Knepley       } else {
1519a5b13a2SMatthew G Knepley         key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0;
1529a5b13a2SMatthew G Knepley         ierr = PetscSortInt(faceSize, (PetscInt *) &key);
1539f074e33SMatthew G Knepley       }
1549a5b13a2SMatthew G Knepley       ierr  = PetscHashIJKLGet(faceTable, key, &f);CHKERRQ(ierr);
1559a5b13a2SMatthew G Knepley       if (f < 0) {
1569a5b13a2SMatthew G Knepley         ierr = PetscHashIJKLAdd(faceTable, key, face);CHKERRQ(ierr);
1579a5b13a2SMatthew G Knepley         f    = face++;
1589a5b13a2SMatthew G Knepley       }
1599a5b13a2SMatthew G Knepley     }
1609a5b13a2SMatthew G Knepley   }
1619a5b13a2SMatthew G Knepley   pEnd[faceDepth] = face;
1629a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr);
1639a5b13a2SMatthew G Knepley   /* Count new points */
1649a5b13a2SMatthew G Knepley   for (d = 0; d <= depth; ++d) {
1659a5b13a2SMatthew G Knepley     numPoints += pEnd[d]-pStart[d];
1669a5b13a2SMatthew G Knepley   }
1679a5b13a2SMatthew G Knepley   ierr = DMPlexSetChart(idm, 0, numPoints);CHKERRQ(ierr);
1689a5b13a2SMatthew G Knepley   /* Set cone sizes */
1699a5b13a2SMatthew G Knepley   for (d = 0; d <= depth; ++d) {
1709a5b13a2SMatthew G Knepley     PetscInt coneSize, p;
1719f074e33SMatthew G Knepley 
1729a5b13a2SMatthew G Knepley     if (d == faceDepth) {
1739a5b13a2SMatthew G Knepley       for (p = pStart[d]; p < pEnd[d]; ++p) {
1749a5b13a2SMatthew G Knepley         /* I see no way to do this if we admit faces of different shapes */
1759a5b13a2SMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr);
1769a5b13a2SMatthew G Knepley       }
177a014044eSMatthew G Knepley     } else if (d == cellDepth) {
178a014044eSMatthew G Knepley       for (p = pStart[d]; p < pEnd[d]; ++p) {
179a014044eSMatthew G Knepley         /* Number of cell faces may be different from number of cell vertices*/
180a014044eSMatthew G Knepley         ierr = DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);CHKERRQ(ierr);
181a014044eSMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr);
182a014044eSMatthew G Knepley       }
1839a5b13a2SMatthew G Knepley     } else {
1849a5b13a2SMatthew G Knepley       for (p = pStart[d]; p < pEnd[d]; ++p) {
1859a5b13a2SMatthew G Knepley         ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
1869a5b13a2SMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr);
1879f074e33SMatthew G Knepley       }
1889f074e33SMatthew G Knepley     }
1899f074e33SMatthew G Knepley   }
1909f074e33SMatthew G Knepley   ierr = DMSetUp(idm);CHKERRQ(ierr);
1919f074e33SMatthew G Knepley   /* Get face cones from subsets of cell vertices */
1929a5b13a2SMatthew 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);
1939a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr);
1949a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLSetMultivalued(faceTable, PETSC_FALSE);CHKERRQ(ierr);
1959a5b13a2SMatthew G Knepley   for (d = depth; d > cellDepth; --d) {
1969a5b13a2SMatthew G Knepley     const PetscInt *cone;
1979a5b13a2SMatthew G Knepley     PetscInt        p;
1989a5b13a2SMatthew G Knepley 
1999a5b13a2SMatthew G Knepley     for (p = pStart[d]; p < pEnd[d]; ++p) {
2009a5b13a2SMatthew G Knepley       ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
2019a5b13a2SMatthew G Knepley       ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr);
2029a5b13a2SMatthew G Knepley       ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr);
2039a5b13a2SMatthew G Knepley       ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr);
2049f074e33SMatthew G Knepley     }
2059a5b13a2SMatthew G Knepley   }
2069a5b13a2SMatthew G Knepley   for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) {
2079f074e33SMatthew G Knepley     const PetscInt *cellFaces;
2089a5b13a2SMatthew G Knepley     PetscInt        numCellFaces, faceSize, cf, f;
2099f074e33SMatthew G Knepley 
2109a5b13a2SMatthew G Knepley     ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
2119a5b13a2SMatthew G Knepley     if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll);
2129f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
2139a5b13a2SMatthew G Knepley       const PetscInt  *cellFace = &cellFaces[cf*faceSize];
2149a5b13a2SMatthew G Knepley       PetscHashIJKLKey key;
2159f074e33SMatthew G Knepley 
2169a5b13a2SMatthew G Knepley       if (faceSize == 2) {
2179a5b13a2SMatthew G Knepley         key.i = PetscMin(cellFace[0], cellFace[1]);
2189a5b13a2SMatthew G Knepley         key.j = PetscMax(cellFace[0], cellFace[1]);
2197fcd4ef0SJed Brown         key.k = 0;
2207fcd4ef0SJed Brown         key.l = 0;
2219a5b13a2SMatthew G Knepley       } else {
2229a5b13a2SMatthew G Knepley         key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0;
2239a5b13a2SMatthew G Knepley         ierr = PetscSortInt(faceSize, (PetscInt *) &key);
2249f074e33SMatthew G Knepley       }
2259a5b13a2SMatthew G Knepley       ierr  = PetscHashIJKLGet(faceTable, key, &f);CHKERRQ(ierr);
2269a5b13a2SMatthew G Knepley       if (f < 0) {
2279a5b13a2SMatthew G Knepley         ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr);
2289a5b13a2SMatthew G Knepley         ierr = PetscHashIJKLAdd(faceTable, key, face);CHKERRQ(ierr);
2299a5b13a2SMatthew G Knepley         f    = face++;
2309f074e33SMatthew G Knepley         ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr);
2319a5b13a2SMatthew G Knepley       } else {
2329a5b13a2SMatthew G Knepley         const PetscInt *cone;
2339a5b13a2SMatthew G Knepley         PetscInt        coneSize, ornt, i, j;
2349f074e33SMatthew G Knepley 
2359a5b13a2SMatthew G Knepley         ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr);
2362e1b13c2SMatthew G. Knepley         /* Orient face: Do not allow reverse orientation at the first vertex */
2379f074e33SMatthew G Knepley         ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr);
2389f074e33SMatthew G Knepley         ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr);
2399a5b13a2SMatthew 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);
2409a5b13a2SMatthew G Knepley         /* - First find the initial vertex */
2419a5b13a2SMatthew G Knepley         for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break;
2429a5b13a2SMatthew G Knepley         /* - Try forward comparison */
2439a5b13a2SMatthew G Knepley         for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break;
2449a5b13a2SMatthew G Knepley         if (j == faceSize) {
2459a5b13a2SMatthew G Knepley           if ((faceSize == 2) && (i == 1)) ornt = -2;
2469a5b13a2SMatthew G Knepley           else                             ornt = i;
2479a5b13a2SMatthew G Knepley         } else {
2489a5b13a2SMatthew G Knepley           /* - Try backward comparison */
2499a5b13a2SMatthew G Knepley           for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break;
2502e1b13c2SMatthew G. Knepley           if (j == faceSize) {
2512e1b13c2SMatthew G. Knepley             if (i == 0) ornt = -faceSize;
2522e1b13c2SMatthew G. Knepley             else        ornt = -(i+1);
2532e1b13c2SMatthew G. Knepley           } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation");
2549a5b13a2SMatthew G Knepley         }
2559a5b13a2SMatthew G Knepley         ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr);
2569f074e33SMatthew G Knepley       }
2579f074e33SMatthew G Knepley     }
2589f074e33SMatthew G Knepley   }
2599a5b13a2SMatthew 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]);
260c907b753SJed Brown   ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr);
2619a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr);
2626551a8c7SMatthew G. Knepley   ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr);
2639a5b13a2SMatthew G Knepley   ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr);
2649a5b13a2SMatthew G Knepley   ierr = DMPlexStratify(idm);CHKERRQ(ierr);
2659f074e33SMatthew G Knepley   PetscFunctionReturn(0);
2669f074e33SMatthew G Knepley }
2679f074e33SMatthew G Knepley 
2689f074e33SMatthew G Knepley #undef __FUNCT__
2699f074e33SMatthew G Knepley #define __FUNCT__ "DMPlexInterpolate"
27080330477SMatthew G. Knepley /*@
27180330477SMatthew G. Knepley   DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc.
27280330477SMatthew G. Knepley 
27380330477SMatthew G. Knepley   Collective on DM
27480330477SMatthew G. Knepley 
27580330477SMatthew G. Knepley   Input Parameter:
276*4e3744c5SMatthew G. Knepley . dm - The DMPlex object with only cells and vertices
27780330477SMatthew G. Knepley 
27880330477SMatthew G. Knepley   Output Parameter:
279*4e3744c5SMatthew G. Knepley . dmInt - The complete DMPlex object
28080330477SMatthew G. Knepley 
28180330477SMatthew G. Knepley   Level: intermediate
28280330477SMatthew G. Knepley 
28380330477SMatthew G. Knepley .keywords: mesh
284*4e3744c5SMatthew G. Knepley .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList()
28580330477SMatthew G. Knepley @*/
2869f074e33SMatthew G Knepley PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt)
2879f074e33SMatthew G Knepley {
2889a5b13a2SMatthew G Knepley   DM             idm, odm = dm;
2892e1b13c2SMatthew G. Knepley   PetscInt       depth, dim, d;
2909f074e33SMatthew G Knepley   PetscErrorCode ierr;
2919f074e33SMatthew G Knepley 
2929f074e33SMatthew G Knepley   PetscFunctionBegin;
2932e1b13c2SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2949f074e33SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
29576b791aaSMatthew G. Knepley   if (dim <= 1) {
29676b791aaSMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
29776b791aaSMatthew G. Knepley     idm  = dm;
29876b791aaSMatthew G. Knepley   }
2999a5b13a2SMatthew G Knepley   for (d = 1; d < dim; ++d) {
3009a5b13a2SMatthew G Knepley     /* Create interpolated mesh */
3019a5b13a2SMatthew G Knepley     ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr);
3029a5b13a2SMatthew G Knepley     ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr);
3039a5b13a2SMatthew G Knepley     ierr = DMPlexSetDimension(idm, dim);CHKERRQ(ierr);
3042e1b13c2SMatthew G. Knepley     if (depth > 0) {ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr);}
3059a5b13a2SMatthew G Knepley     if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);}
3069a5b13a2SMatthew G Knepley     odm  = idm;
3079f074e33SMatthew G Knepley   }
3089a5b13a2SMatthew G Knepley   *dmInt = idm;
3099f074e33SMatthew G Knepley   PetscFunctionReturn(0);
3109f074e33SMatthew G Knepley }
31107b0a1fcSMatthew G Knepley 
31207b0a1fcSMatthew G Knepley #undef __FUNCT__
31307b0a1fcSMatthew G Knepley #define __FUNCT__ "DMPlexCopyCoordinates"
31480330477SMatthew G. Knepley /*@
31580330477SMatthew G. Knepley   DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices
31680330477SMatthew G. Knepley 
31780330477SMatthew G. Knepley   Collective on DM
31880330477SMatthew G. Knepley 
31980330477SMatthew G. Knepley   Input Parameter:
32080330477SMatthew G. Knepley . dmA - The DMPlex object with initial coordinates
32180330477SMatthew G. Knepley 
32280330477SMatthew G. Knepley   Output Parameter:
32380330477SMatthew G. Knepley . dmB - The DMPlex object with copied coordinates
32480330477SMatthew G. Knepley 
32580330477SMatthew G. Knepley   Level: intermediate
32680330477SMatthew G. Knepley 
32780330477SMatthew G. Knepley   Note: This is typically used when adding pieces other than vertices to a mesh
32880330477SMatthew G. Knepley 
32980330477SMatthew G. Knepley .keywords: mesh
33080330477SMatthew G. Knepley .seealso: DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMPlexGetCoordinateSection()
33180330477SMatthew G. Knepley @*/
33207b0a1fcSMatthew G Knepley PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB)
33307b0a1fcSMatthew G Knepley {
33407b0a1fcSMatthew G Knepley   Vec            coordinatesA, coordinatesB;
33507b0a1fcSMatthew G Knepley   PetscSection   coordSectionA, coordSectionB;
33607b0a1fcSMatthew G Knepley   PetscScalar   *coordsA, *coordsB;
33707b0a1fcSMatthew G Knepley   PetscInt       spaceDim, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d;
33807b0a1fcSMatthew G Knepley   PetscErrorCode ierr;
33907b0a1fcSMatthew G Knepley 
34007b0a1fcSMatthew G Knepley   PetscFunctionBegin;
34176b791aaSMatthew G. Knepley   if (dmA == dmB) PetscFunctionReturn(0);
34207b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr);
34307b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr);
34407b0a1fcSMatthew 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);
34507b0a1fcSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr);
34607b0a1fcSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr);
34707b0a1fcSMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr);
34807b0a1fcSMatthew G Knepley   ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr);
34907b0a1fcSMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr);
35007b0a1fcSMatthew G Knepley   ierr = PetscSectionSetChart(coordSectionB, vStartB, vEndB);CHKERRQ(ierr);
35107b0a1fcSMatthew G Knepley   for (v = vStartB; v < vEndB; ++v) {
35207b0a1fcSMatthew G Knepley     ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr);
35307b0a1fcSMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr);
35407b0a1fcSMatthew G Knepley   }
35507b0a1fcSMatthew G Knepley   ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr);
35607b0a1fcSMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr);
35707b0a1fcSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr);
35807b0a1fcSMatthew G Knepley   ierr = VecCreate(PetscObjectComm((PetscObject) dmB), &coordinatesB);CHKERRQ(ierr);
35907b0a1fcSMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr);
36007b0a1fcSMatthew G Knepley   ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr);
36107b0a1fcSMatthew G Knepley   ierr = VecSetFromOptions(coordinatesB);CHKERRQ(ierr);
36207b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr);
36307b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr);
36407b0a1fcSMatthew G Knepley   for (v = 0; v < vEndB-vStartB; ++v) {
36507b0a1fcSMatthew G Knepley     for (d = 0; d < spaceDim; ++d) {
36607b0a1fcSMatthew G Knepley       coordsB[v*spaceDim+d] = coordsA[v*spaceDim+d];
36707b0a1fcSMatthew G Knepley     }
36807b0a1fcSMatthew G Knepley   }
36907b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr);
37007b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr);
37107b0a1fcSMatthew G Knepley   ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr);
37207b0a1fcSMatthew G Knepley   ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr);
37307b0a1fcSMatthew G Knepley   PetscFunctionReturn(0);
37407b0a1fcSMatthew G Knepley }
375*4e3744c5SMatthew G. Knepley 
376*4e3744c5SMatthew G. Knepley #undef __FUNCT__
377*4e3744c5SMatthew G. Knepley #define __FUNCT__ "DMPlexUninterpolate"
378*4e3744c5SMatthew G. Knepley /*@
379*4e3744c5SMatthew G. Knepley   DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh
380*4e3744c5SMatthew G. Knepley 
381*4e3744c5SMatthew G. Knepley   Collective on DM
382*4e3744c5SMatthew G. Knepley 
383*4e3744c5SMatthew G. Knepley   Input Parameter:
384*4e3744c5SMatthew G. Knepley . dm - The complete DMPlex object
385*4e3744c5SMatthew G. Knepley 
386*4e3744c5SMatthew G. Knepley   Output Parameter:
387*4e3744c5SMatthew G. Knepley . dmUnint - The DMPlex object with only cells and vertices
388*4e3744c5SMatthew G. Knepley 
389*4e3744c5SMatthew G. Knepley   Level: intermediate
390*4e3744c5SMatthew G. Knepley 
391*4e3744c5SMatthew G. Knepley .keywords: mesh
392*4e3744c5SMatthew G. Knepley .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList()
393*4e3744c5SMatthew G. Knepley @*/
394*4e3744c5SMatthew G. Knepley PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint)
395*4e3744c5SMatthew G. Knepley {
396*4e3744c5SMatthew G. Knepley   DM             udm;
397*4e3744c5SMatthew G. Knepley   PetscInt       dim, vStart, vEnd, cStart, cEnd, c, maxConeSize, *cone;
398*4e3744c5SMatthew G. Knepley   PetscErrorCode ierr;
399*4e3744c5SMatthew G. Knepley 
400*4e3744c5SMatthew G. Knepley   PetscFunctionBegin;
401*4e3744c5SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
402*4e3744c5SMatthew G. Knepley   if (dim <= 1) {
403*4e3744c5SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
404*4e3744c5SMatthew G. Knepley     udm  = dm;
405*4e3744c5SMatthew G. Knepley   }
406*4e3744c5SMatthew G. Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
407*4e3744c5SMatthew G. Knepley   ierr = DMPlexGetDepthStratum(udm, 0, &vStart, &vEnd);CHKERRQ(ierr);
408*4e3744c5SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(udm, 0, &cStart, &cEnd);CHKERRQ(ierr);
409*4e3744c5SMatthew G. Knepley   ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr);
410*4e3744c5SMatthew G. Knepley   ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr);
411*4e3744c5SMatthew G. Knepley   ierr = DMPlexSetDimension(udm, dim);CHKERRQ(ierr);
412*4e3744c5SMatthew G. Knepley   ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr);
413*4e3744c5SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
414*4e3744c5SMatthew G. Knepley     PetscInt *closure, closureSize, cl, coneSize = 0;
415*4e3744c5SMatthew G. Knepley 
416*4e3744c5SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
417*4e3744c5SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
418*4e3744c5SMatthew G. Knepley       const PetscInt p = closure[cl];
419*4e3744c5SMatthew G. Knepley 
420*4e3744c5SMatthew G. Knepley       if ((p >= vStart) && (p < vEnd)) ++coneSize;
421*4e3744c5SMatthew G. Knepley     }
422*4e3744c5SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
423*4e3744c5SMatthew G. Knepley     ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr);
424*4e3744c5SMatthew G. Knepley   }
425*4e3744c5SMatthew G. Knepley   ierr = DMSetUp(udm);CHKERRQ(ierr);
426*4e3744c5SMatthew G. Knepley   ierr = PetscMalloc(maxConeSize * sizeof(PetscInt), &cone);CHKERRQ(ierr);
427*4e3744c5SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
428*4e3744c5SMatthew G. Knepley     PetscInt *closure, closureSize, cl, coneSize = 0;
429*4e3744c5SMatthew G. Knepley 
430*4e3744c5SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
431*4e3744c5SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
432*4e3744c5SMatthew G. Knepley       const PetscInt p = closure[cl];
433*4e3744c5SMatthew G. Knepley 
434*4e3744c5SMatthew G. Knepley       if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p;
435*4e3744c5SMatthew G. Knepley     }
436*4e3744c5SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
437*4e3744c5SMatthew G. Knepley     ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr);
438*4e3744c5SMatthew G. Knepley   }
439*4e3744c5SMatthew G. Knepley   ierr = PetscFree(cone);CHKERRQ(ierr);
440*4e3744c5SMatthew G. Knepley   ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr);
441*4e3744c5SMatthew G. Knepley   ierr = DMPlexStratify(udm);CHKERRQ(ierr);
442*4e3744c5SMatthew G. Knepley   *dmUnint = udm;
443*4e3744c5SMatthew G. Knepley   PetscFunctionReturn(0);
444*4e3744c5SMatthew G. Knepley }
445