xref: /petsc/src/dm/impls/plex/plexinterpolate.c (revision 439ece16d4e67c1951eb9b293e2dfbd1258932d1)
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 */
9*439ece16SMatthew 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       *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);
22*439ece16SMatthew G. Knepley   ierr = DMPlexGetRawFaces_Internal(dm, dim, coneSize, cone, numFaces, faceSize, faces);CHKERRQ(ierr);
23*439ece16SMatthew G. Knepley   PetscFunctionReturn(0);
24*439ece16SMatthew G. Knepley }
25*439ece16SMatthew G. Knepley 
26*439ece16SMatthew G. Knepley #undef __FUNCT__
27*439ece16SMatthew G. Knepley #define __FUNCT__ "DMPlexRestoreFaces_Internal"
28*439ece16SMatthew G. Knepley /*
29*439ece16SMatthew G. Knepley   DMPlexRestoreFaces_Internal - Restores the array
30*439ece16SMatthew G. Knepley */
31*439ece16SMatthew G. Knepley PetscErrorCode DMPlexRestoreFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[])
32*439ece16SMatthew G. Knepley {
33*439ece16SMatthew G. Knepley   PetscErrorCode  ierr;
34*439ece16SMatthew G. Knepley 
35*439ece16SMatthew G. Knepley   PetscFunctionBegin;
36*439ece16SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, 0, PETSC_INT, faces);CHKERRQ(ierr);
37*439ece16SMatthew G. Knepley   PetscFunctionReturn(0);
38*439ece16SMatthew G. Knepley }
39*439ece16SMatthew G. Knepley 
40*439ece16SMatthew G. Knepley #undef __FUNCT__
41*439ece16SMatthew G. Knepley #define __FUNCT__ "DMPlexGetRawFaces_Internal"
42*439ece16SMatthew G. Knepley /*
43*439ece16SMatthew G. Knepley   DMPlexGetRawFaces_Internal - Gets groups of vertices that correspond to faces for the given cone
44*439ece16SMatthew G. Knepley */
45*439ece16SMatthew G. Knepley PetscErrorCode DMPlexGetRawFaces_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[])
46*439ece16SMatthew G. Knepley {
47*439ece16SMatthew G. Knepley   PetscInt       *facesTmp;
48*439ece16SMatthew G. Knepley   PetscInt        maxConeSize, maxSupportSize;
49*439ece16SMatthew G. Knepley   PetscErrorCode  ierr;
50*439ece16SMatthew G. Knepley 
51*439ece16SMatthew G. Knepley   PetscFunctionBegin;
52*439ece16SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
53*439ece16SMatthew G. Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
54*439ece16SMatthew G. Knepley   ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), PETSC_INT, &facesTmp);CHKERRQ(ierr);
559f074e33SMatthew G Knepley   switch (dim) {
569f074e33SMatthew G Knepley   case 2:
579f074e33SMatthew G Knepley     switch (coneSize) {
589f074e33SMatthew G Knepley     case 3:
599a5b13a2SMatthew G Knepley       if (faces) {
609a5b13a2SMatthew G Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
619a5b13a2SMatthew G Knepley         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
629a5b13a2SMatthew G Knepley         facesTmp[4] = cone[2]; facesTmp[5] = cone[0];
639a5b13a2SMatthew G Knepley         *faces = facesTmp;
649a5b13a2SMatthew G Knepley       }
659a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 3;
669a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 2;
679f074e33SMatthew G Knepley       break;
689f074e33SMatthew G Knepley     case 4:
699a5b13a2SMatthew G Knepley       /* Vertices follow right hand rule */
709a5b13a2SMatthew G Knepley       if (faces) {
719a5b13a2SMatthew G Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
729a5b13a2SMatthew G Knepley         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
739a5b13a2SMatthew G Knepley         facesTmp[4] = cone[2]; facesTmp[5] = cone[3];
749a5b13a2SMatthew G Knepley         facesTmp[6] = cone[3]; facesTmp[7] = cone[0];
759a5b13a2SMatthew G Knepley         *faces = facesTmp;
769a5b13a2SMatthew G Knepley       }
779a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 4;
789a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 2;
799a5b13a2SMatthew G Knepley       if (faces)    *faces            = facesTmp;
809f074e33SMatthew G Knepley       break;
819f074e33SMatthew G Knepley     default:
829f074e33SMatthew G Knepley       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
839f074e33SMatthew G Knepley     }
849f074e33SMatthew G Knepley     break;
859f074e33SMatthew G Knepley   case 3:
869f074e33SMatthew G Knepley     switch (coneSize) {
879f074e33SMatthew G Knepley     case 3:
889a5b13a2SMatthew G Knepley       if (faces) {
899a5b13a2SMatthew G Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
909a5b13a2SMatthew G Knepley         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
919a5b13a2SMatthew G Knepley         facesTmp[4] = cone[2]; facesTmp[5] = cone[0];
929a5b13a2SMatthew G Knepley         *faces = facesTmp;
939a5b13a2SMatthew G Knepley       }
949a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 3;
959a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 2;
969a5b13a2SMatthew G Knepley       if (faces)    *faces            = facesTmp;
979f074e33SMatthew G Knepley       break;
989f074e33SMatthew G Knepley     case 4:
992e1b13c2SMatthew G. Knepley       /* Vertices of first face follow right hand rule and normal points away from last vertex */
1009a5b13a2SMatthew G Knepley       if (faces) {
1012e1b13c2SMatthew G. Knepley         facesTmp[0] = cone[0]; facesTmp[1]  = cone[1]; facesTmp[2]  = cone[2];
1022e1b13c2SMatthew G. Knepley         facesTmp[3] = cone[0]; facesTmp[4]  = cone[3]; facesTmp[5]  = cone[1];
1032e1b13c2SMatthew G. Knepley         facesTmp[6] = cone[0]; facesTmp[7]  = cone[2]; facesTmp[8]  = cone[3];
1042e1b13c2SMatthew G. Knepley         facesTmp[9] = cone[2]; facesTmp[10] = cone[1]; facesTmp[11] = cone[3];
1059a5b13a2SMatthew G Knepley         *faces = facesTmp;
1069a5b13a2SMatthew G Knepley       }
1079a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 4;
1089a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 3;
1099a5b13a2SMatthew G Knepley       if (faces)    *faces            = facesTmp;
1109a5b13a2SMatthew G Knepley       break;
1119a5b13a2SMatthew G Knepley     case 8:
1129a5b13a2SMatthew G Knepley       if (faces) {
1132e1b13c2SMatthew G. Knepley         facesTmp[0]  = cone[0]; facesTmp[1]  = cone[1]; facesTmp[2]  = cone[2]; facesTmp[3]  = cone[3];
114a014044eSMatthew G Knepley         facesTmp[4]  = cone[4]; facesTmp[5]  = cone[5]; facesTmp[6]  = cone[6]; facesTmp[7]  = cone[7];
1152e1b13c2SMatthew G. Knepley         facesTmp[8]  = cone[0]; facesTmp[9]  = cone[3]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4];
1162e1b13c2SMatthew G. Knepley         facesTmp[12] = cone[2]; facesTmp[13] = cone[1]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6];
1172e1b13c2SMatthew G. Knepley         facesTmp[16] = cone[3]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5];
1182e1b13c2SMatthew G. Knepley         facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[1];
1199a5b13a2SMatthew G Knepley         *faces = facesTmp;
1209a5b13a2SMatthew G Knepley       }
1219a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 6;
1229a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 4;
1239f074e33SMatthew G Knepley       break;
1249f074e33SMatthew G Knepley     default:
1259f074e33SMatthew G Knepley       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
1269f074e33SMatthew G Knepley     }
1279f074e33SMatthew G Knepley     break;
1289f074e33SMatthew G Knepley   default:
1299f074e33SMatthew G Knepley     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim);
1309f074e33SMatthew G Knepley   }
1319f074e33SMatthew G Knepley   PetscFunctionReturn(0);
1329f074e33SMatthew G Knepley }
1339f074e33SMatthew G Knepley 
1349f074e33SMatthew G Knepley #undef __FUNCT__
1359a5b13a2SMatthew G Knepley #define __FUNCT__ "DMPlexInterpolateFaces_Internal"
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;
1469a5b13a2SMatthew G Knepley   ierr = DMPlexGetDimension(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;
1549a5b13a2SMatthew G Knepley   ierr = PetscMalloc2(depth+1,PetscInt,&pStart,depth+1,PetscInt,&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   ierr = PetscHashIJKLSetMultivalued(faceTable, PETSC_FALSE);CHKERRQ(ierr);
1679a5b13a2SMatthew G Knepley   for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) {
1689f074e33SMatthew G Knepley     const PetscInt *cellFaces;
1699a5b13a2SMatthew G Knepley     PetscInt        numCellFaces, faceSize, cf, f;
1709f074e33SMatthew G Knepley 
1719a5b13a2SMatthew G Knepley     ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
1729a5b13a2SMatthew G Knepley     if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll);
1739f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
1749a5b13a2SMatthew G Knepley       const PetscInt  *cellFace = &cellFaces[cf*faceSize];
1759a5b13a2SMatthew G Knepley       PetscHashIJKLKey key;
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;
1849a5b13a2SMatthew G Knepley         ierr = PetscSortInt(faceSize, (PetscInt *) &key);
1859f074e33SMatthew G Knepley       }
1869a5b13a2SMatthew G Knepley       ierr  = PetscHashIJKLGet(faceTable, key, &f);CHKERRQ(ierr);
1879a5b13a2SMatthew G Knepley       if (f < 0) {
1889a5b13a2SMatthew G Knepley         ierr = PetscHashIJKLAdd(faceTable, key, face);CHKERRQ(ierr);
1899a5b13a2SMatthew G Knepley         f    = face++;
1909a5b13a2SMatthew G Knepley       }
1919a5b13a2SMatthew G Knepley     }
192*439ece16SMatthew G. Knepley     ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
1939a5b13a2SMatthew G Knepley   }
1949a5b13a2SMatthew G Knepley   pEnd[faceDepth] = face;
1959a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr);
1969a5b13a2SMatthew G Knepley   /* Count new points */
1979a5b13a2SMatthew G Knepley   for (d = 0; d <= depth; ++d) {
1989a5b13a2SMatthew G Knepley     numPoints += pEnd[d]-pStart[d];
1999a5b13a2SMatthew G Knepley   }
2009a5b13a2SMatthew G Knepley   ierr = DMPlexSetChart(idm, 0, numPoints);CHKERRQ(ierr);
2019a5b13a2SMatthew G Knepley   /* Set cone sizes */
2029a5b13a2SMatthew G Knepley   for (d = 0; d <= depth; ++d) {
2039a5b13a2SMatthew G Knepley     PetscInt coneSize, p;
2049f074e33SMatthew G Knepley 
2059a5b13a2SMatthew G Knepley     if (d == faceDepth) {
2069a5b13a2SMatthew G Knepley       for (p = pStart[d]; p < pEnd[d]; ++p) {
2079a5b13a2SMatthew G Knepley         /* I see no way to do this if we admit faces of different shapes */
2089a5b13a2SMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr);
2099a5b13a2SMatthew G Knepley       }
210a014044eSMatthew G Knepley     } else if (d == cellDepth) {
211a014044eSMatthew G Knepley       for (p = pStart[d]; p < pEnd[d]; ++p) {
212a014044eSMatthew G Knepley         /* Number of cell faces may be different from number of cell vertices*/
213a014044eSMatthew G Knepley         ierr = DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);CHKERRQ(ierr);
214a014044eSMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr);
215a014044eSMatthew G Knepley       }
2169a5b13a2SMatthew G Knepley     } else {
2179a5b13a2SMatthew G Knepley       for (p = pStart[d]; p < pEnd[d]; ++p) {
2189a5b13a2SMatthew G Knepley         ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
2199a5b13a2SMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr);
2209f074e33SMatthew G Knepley       }
2219f074e33SMatthew G Knepley     }
2229f074e33SMatthew G Knepley   }
2239f074e33SMatthew G Knepley   ierr = DMSetUp(idm);CHKERRQ(ierr);
2249f074e33SMatthew G Knepley   /* Get face cones from subsets of cell vertices */
2259a5b13a2SMatthew 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);
2269a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr);
2279a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLSetMultivalued(faceTable, PETSC_FALSE);CHKERRQ(ierr);
2289a5b13a2SMatthew G Knepley   for (d = depth; d > cellDepth; --d) {
2299a5b13a2SMatthew G Knepley     const PetscInt *cone;
2309a5b13a2SMatthew G Knepley     PetscInt        p;
2319a5b13a2SMatthew G Knepley 
2329a5b13a2SMatthew G Knepley     for (p = pStart[d]; p < pEnd[d]; ++p) {
2339a5b13a2SMatthew G Knepley       ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
2349a5b13a2SMatthew G Knepley       ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr);
2359a5b13a2SMatthew G Knepley       ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr);
2369a5b13a2SMatthew G Knepley       ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr);
2379f074e33SMatthew G Knepley     }
2389a5b13a2SMatthew G Knepley   }
2399a5b13a2SMatthew G Knepley   for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) {
2409f074e33SMatthew G Knepley     const PetscInt *cellFaces;
2419a5b13a2SMatthew G Knepley     PetscInt        numCellFaces, faceSize, cf, f;
2429f074e33SMatthew G Knepley 
2439a5b13a2SMatthew G Knepley     ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
2449a5b13a2SMatthew G Knepley     if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll);
2459f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
2469a5b13a2SMatthew G Knepley       const PetscInt  *cellFace = &cellFaces[cf*faceSize];
2479a5b13a2SMatthew G Knepley       PetscHashIJKLKey key;
2489f074e33SMatthew G Knepley 
2499a5b13a2SMatthew G Knepley       if (faceSize == 2) {
2509a5b13a2SMatthew G Knepley         key.i = PetscMin(cellFace[0], cellFace[1]);
2519a5b13a2SMatthew G Knepley         key.j = PetscMax(cellFace[0], cellFace[1]);
2527fcd4ef0SJed Brown         key.k = 0;
2537fcd4ef0SJed Brown         key.l = 0;
2549a5b13a2SMatthew G Knepley       } else {
2559a5b13a2SMatthew G Knepley         key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0;
2569a5b13a2SMatthew G Knepley         ierr = PetscSortInt(faceSize, (PetscInt *) &key);
2579f074e33SMatthew G Knepley       }
2589a5b13a2SMatthew G Knepley       ierr  = PetscHashIJKLGet(faceTable, key, &f);CHKERRQ(ierr);
2599a5b13a2SMatthew G Knepley       if (f < 0) {
2609a5b13a2SMatthew G Knepley         ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr);
2619a5b13a2SMatthew G Knepley         ierr = PetscHashIJKLAdd(faceTable, key, face);CHKERRQ(ierr);
2629a5b13a2SMatthew G Knepley         f    = face++;
2639f074e33SMatthew G Knepley         ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr);
2649a5b13a2SMatthew G Knepley       } else {
2659a5b13a2SMatthew G Knepley         const PetscInt *cone;
2669a5b13a2SMatthew G Knepley         PetscInt        coneSize, ornt, i, j;
2679f074e33SMatthew G Knepley 
2689a5b13a2SMatthew G Knepley         ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr);
2692e1b13c2SMatthew G. Knepley         /* Orient face: Do not allow reverse orientation at the first vertex */
2709f074e33SMatthew G Knepley         ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr);
2719f074e33SMatthew G Knepley         ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr);
2729a5b13a2SMatthew 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);
2739a5b13a2SMatthew G Knepley         /* - First find the initial vertex */
2749a5b13a2SMatthew G Knepley         for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break;
2759a5b13a2SMatthew G Knepley         /* - Try forward comparison */
2769a5b13a2SMatthew G Knepley         for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break;
2779a5b13a2SMatthew G Knepley         if (j == faceSize) {
2789a5b13a2SMatthew G Knepley           if ((faceSize == 2) && (i == 1)) ornt = -2;
2799a5b13a2SMatthew G Knepley           else                             ornt = i;
2809a5b13a2SMatthew G Knepley         } else {
2819a5b13a2SMatthew G Knepley           /* - Try backward comparison */
2829a5b13a2SMatthew G Knepley           for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break;
2832e1b13c2SMatthew G. Knepley           if (j == faceSize) {
2842e1b13c2SMatthew G. Knepley             if (i == 0) ornt = -faceSize;
2852e1b13c2SMatthew G. Knepley             else        ornt = -(i+1);
2862e1b13c2SMatthew G. Knepley           } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation");
2879a5b13a2SMatthew G Knepley         }
2889a5b13a2SMatthew G Knepley         ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr);
2899f074e33SMatthew G Knepley       }
2909f074e33SMatthew G Knepley     }
291*439ece16SMatthew G. Knepley     ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
2929f074e33SMatthew G Knepley   }
2939a5b13a2SMatthew 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]);
294c907b753SJed Brown   ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr);
2959a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr);
2966551a8c7SMatthew G. Knepley   ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr);
2979a5b13a2SMatthew G Knepley   ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr);
2989a5b13a2SMatthew G Knepley   ierr = DMPlexStratify(idm);CHKERRQ(ierr);
2999f074e33SMatthew G Knepley   PetscFunctionReturn(0);
3009f074e33SMatthew G Knepley }
3019f074e33SMatthew G Knepley 
3029f074e33SMatthew G Knepley #undef __FUNCT__
3039f074e33SMatthew G Knepley #define __FUNCT__ "DMPlexInterpolate"
30480330477SMatthew G. Knepley /*@
30580330477SMatthew G. Knepley   DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc.
30680330477SMatthew G. Knepley 
30780330477SMatthew G. Knepley   Collective on DM
30880330477SMatthew G. Knepley 
30980330477SMatthew G. Knepley   Input Parameter:
3104e3744c5SMatthew G. Knepley . dm - The DMPlex object with only cells and vertices
31180330477SMatthew G. Knepley 
31280330477SMatthew G. Knepley   Output Parameter:
3134e3744c5SMatthew G. Knepley . dmInt - The complete DMPlex object
31480330477SMatthew G. Knepley 
31580330477SMatthew G. Knepley   Level: intermediate
31680330477SMatthew G. Knepley 
31780330477SMatthew G. Knepley .keywords: mesh
3184e3744c5SMatthew G. Knepley .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList()
31980330477SMatthew G. Knepley @*/
3209f074e33SMatthew G Knepley PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt)
3219f074e33SMatthew G Knepley {
3229a5b13a2SMatthew G Knepley   DM             idm, odm = dm;
3232e1b13c2SMatthew G. Knepley   PetscInt       depth, dim, d;
3249f074e33SMatthew G Knepley   PetscErrorCode ierr;
3259f074e33SMatthew G Knepley 
3269f074e33SMatthew G Knepley   PetscFunctionBegin;
3272e1b13c2SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
3289f074e33SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
32976b791aaSMatthew G. Knepley   if (dim <= 1) {
33076b791aaSMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
33176b791aaSMatthew G. Knepley     idm  = dm;
33276b791aaSMatthew G. Knepley   }
3339a5b13a2SMatthew G Knepley   for (d = 1; d < dim; ++d) {
3349a5b13a2SMatthew G Knepley     /* Create interpolated mesh */
3359a5b13a2SMatthew G Knepley     ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr);
3369a5b13a2SMatthew G Knepley     ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr);
3379a5b13a2SMatthew G Knepley     ierr = DMPlexSetDimension(idm, dim);CHKERRQ(ierr);
3382e1b13c2SMatthew G. Knepley     if (depth > 0) {ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr);}
3399a5b13a2SMatthew G Knepley     if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);}
3409a5b13a2SMatthew G Knepley     odm  = idm;
3419f074e33SMatthew G Knepley   }
3429a5b13a2SMatthew G Knepley   *dmInt = idm;
3439f074e33SMatthew G Knepley   PetscFunctionReturn(0);
3449f074e33SMatthew G Knepley }
34507b0a1fcSMatthew G Knepley 
34607b0a1fcSMatthew G Knepley #undef __FUNCT__
34707b0a1fcSMatthew G Knepley #define __FUNCT__ "DMPlexCopyCoordinates"
34880330477SMatthew G. Knepley /*@
34980330477SMatthew G. Knepley   DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices
35080330477SMatthew G. Knepley 
35180330477SMatthew G. Knepley   Collective on DM
35280330477SMatthew G. Knepley 
35380330477SMatthew G. Knepley   Input Parameter:
35480330477SMatthew G. Knepley . dmA - The DMPlex object with initial coordinates
35580330477SMatthew G. Knepley 
35680330477SMatthew G. Knepley   Output Parameter:
35780330477SMatthew G. Knepley . dmB - The DMPlex object with copied coordinates
35880330477SMatthew G. Knepley 
35980330477SMatthew G. Knepley   Level: intermediate
36080330477SMatthew G. Knepley 
36180330477SMatthew G. Knepley   Note: This is typically used when adding pieces other than vertices to a mesh
36280330477SMatthew G. Knepley 
36380330477SMatthew G. Knepley .keywords: mesh
36480330477SMatthew G. Knepley .seealso: DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMPlexGetCoordinateSection()
36580330477SMatthew G. Knepley @*/
36607b0a1fcSMatthew G Knepley PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB)
36707b0a1fcSMatthew G Knepley {
36807b0a1fcSMatthew G Knepley   Vec            coordinatesA, coordinatesB;
36907b0a1fcSMatthew G Knepley   PetscSection   coordSectionA, coordSectionB;
37007b0a1fcSMatthew G Knepley   PetscScalar   *coordsA, *coordsB;
37107b0a1fcSMatthew G Knepley   PetscInt       spaceDim, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d;
37207b0a1fcSMatthew G Knepley   PetscErrorCode ierr;
37307b0a1fcSMatthew G Knepley 
37407b0a1fcSMatthew G Knepley   PetscFunctionBegin;
37576b791aaSMatthew G. Knepley   if (dmA == dmB) PetscFunctionReturn(0);
37607b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr);
37707b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr);
37807b0a1fcSMatthew 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);
37907b0a1fcSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr);
38007b0a1fcSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr);
38107b0a1fcSMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr);
38207b0a1fcSMatthew G Knepley   ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr);
38307b0a1fcSMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr);
38407b0a1fcSMatthew G Knepley   ierr = PetscSectionSetChart(coordSectionB, vStartB, vEndB);CHKERRQ(ierr);
38507b0a1fcSMatthew G Knepley   for (v = vStartB; v < vEndB; ++v) {
38607b0a1fcSMatthew G Knepley     ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr);
38707b0a1fcSMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr);
38807b0a1fcSMatthew G Knepley   }
38907b0a1fcSMatthew G Knepley   ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr);
39007b0a1fcSMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr);
39107b0a1fcSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr);
39207b0a1fcSMatthew G Knepley   ierr = VecCreate(PetscObjectComm((PetscObject) dmB), &coordinatesB);CHKERRQ(ierr);
39307b0a1fcSMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr);
39407b0a1fcSMatthew G Knepley   ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr);
39507b0a1fcSMatthew G Knepley   ierr = VecSetFromOptions(coordinatesB);CHKERRQ(ierr);
39607b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr);
39707b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr);
39807b0a1fcSMatthew G Knepley   for (v = 0; v < vEndB-vStartB; ++v) {
39907b0a1fcSMatthew G Knepley     for (d = 0; d < spaceDim; ++d) {
40007b0a1fcSMatthew G Knepley       coordsB[v*spaceDim+d] = coordsA[v*spaceDim+d];
40107b0a1fcSMatthew G Knepley     }
40207b0a1fcSMatthew G Knepley   }
40307b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr);
40407b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr);
40507b0a1fcSMatthew G Knepley   ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr);
40607b0a1fcSMatthew G Knepley   ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr);
40707b0a1fcSMatthew G Knepley   PetscFunctionReturn(0);
40807b0a1fcSMatthew G Knepley }
4094e3744c5SMatthew G. Knepley 
4104e3744c5SMatthew G. Knepley #undef __FUNCT__
4114e3744c5SMatthew G. Knepley #define __FUNCT__ "DMPlexUninterpolate"
4124e3744c5SMatthew G. Knepley /*@
4134e3744c5SMatthew G. Knepley   DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh
4144e3744c5SMatthew G. Knepley 
4154e3744c5SMatthew G. Knepley   Collective on DM
4164e3744c5SMatthew G. Knepley 
4174e3744c5SMatthew G. Knepley   Input Parameter:
4184e3744c5SMatthew G. Knepley . dm - The complete DMPlex object
4194e3744c5SMatthew G. Knepley 
4204e3744c5SMatthew G. Knepley   Output Parameter:
4214e3744c5SMatthew G. Knepley . dmUnint - The DMPlex object with only cells and vertices
4224e3744c5SMatthew G. Knepley 
4234e3744c5SMatthew G. Knepley   Level: intermediate
4244e3744c5SMatthew G. Knepley 
4254e3744c5SMatthew G. Knepley .keywords: mesh
4264e3744c5SMatthew G. Knepley .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList()
4274e3744c5SMatthew G. Knepley @*/
4284e3744c5SMatthew G. Knepley PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint)
4294e3744c5SMatthew G. Knepley {
4304e3744c5SMatthew G. Knepley   DM             udm;
431595d4abbSMatthew G. Knepley   PetscInt       dim, vStart, vEnd, cStart, cEnd, c, maxConeSize = 0, *cone;
4324e3744c5SMatthew G. Knepley   PetscErrorCode ierr;
4334e3744c5SMatthew G. Knepley 
4344e3744c5SMatthew G. Knepley   PetscFunctionBegin;
4354e3744c5SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
4364e3744c5SMatthew G. Knepley   if (dim <= 1) {
4374e3744c5SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
438595d4abbSMatthew G. Knepley     *dmUnint = dm;
439595d4abbSMatthew G. Knepley     PetscFunctionReturn(0);
4404e3744c5SMatthew G. Knepley   }
441595d4abbSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
442595d4abbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
4434e3744c5SMatthew G. Knepley   ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr);
4444e3744c5SMatthew G. Knepley   ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr);
4454e3744c5SMatthew G. Knepley   ierr = DMPlexSetDimension(udm, dim);CHKERRQ(ierr);
4464e3744c5SMatthew G. Knepley   ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr);
4474e3744c5SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
448595d4abbSMatthew G. Knepley     PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
4494e3744c5SMatthew G. Knepley 
4504e3744c5SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
4514e3744c5SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
4524e3744c5SMatthew G. Knepley       const PetscInt p = closure[cl];
4534e3744c5SMatthew G. Knepley 
4544e3744c5SMatthew G. Knepley       if ((p >= vStart) && (p < vEnd)) ++coneSize;
4554e3744c5SMatthew G. Knepley     }
4564e3744c5SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
4574e3744c5SMatthew G. Knepley     ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr);
458595d4abbSMatthew G. Knepley     maxConeSize = PetscMax(maxConeSize, coneSize);
4594e3744c5SMatthew G. Knepley   }
4604e3744c5SMatthew G. Knepley   ierr = DMSetUp(udm);CHKERRQ(ierr);
4614e3744c5SMatthew G. Knepley   ierr = PetscMalloc(maxConeSize * sizeof(PetscInt), &cone);CHKERRQ(ierr);
4624e3744c5SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
463595d4abbSMatthew G. Knepley     PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
4644e3744c5SMatthew G. Knepley 
4654e3744c5SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
4664e3744c5SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
4674e3744c5SMatthew G. Knepley       const PetscInt p = closure[cl];
4684e3744c5SMatthew G. Knepley 
4694e3744c5SMatthew G. Knepley       if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p;
4704e3744c5SMatthew G. Knepley     }
4714e3744c5SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
4724e3744c5SMatthew G. Knepley     ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr);
4734e3744c5SMatthew G. Knepley   }
4744e3744c5SMatthew G. Knepley   ierr = PetscFree(cone);CHKERRQ(ierr);
4754e3744c5SMatthew G. Knepley   ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr);
4764e3744c5SMatthew G. Knepley   ierr = DMPlexStratify(udm);CHKERRQ(ierr);
4774e3744c5SMatthew G. Knepley   *dmUnint = udm;
4784e3744c5SMatthew G. Knepley   PetscFunctionReturn(0);
4794e3744c5SMatthew G. Knepley }
480