xref: /petsc/src/dm/impls/plex/plexinterpolate.c (revision 9a5b13a2b83addd8b8a14a54302678debae28b96)
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__
5*9a5b13a2SMatthew G Knepley #define __FUNCT__ "DMPlexGetFaces_Internal"
69f074e33SMatthew G Knepley /*
7*9a5b13a2SMatthew G Knepley   DMPlexGetFaces_Internal - Gets groups of vertices that correspond to faces for the given cell
89f074e33SMatthew G Knepley */
9*9a5b13a2SMatthew 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;
12*9a5b13a2SMatthew G Knepley   PetscInt       *facesTmp;
13*9a5b13a2SMatthew 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);
18*9a5b13a2SMatthew G Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
19*9a5b13a2SMatthew 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:
26*9a5b13a2SMatthew G Knepley       if (faces) {
27*9a5b13a2SMatthew G Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
28*9a5b13a2SMatthew G Knepley         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
29*9a5b13a2SMatthew G Knepley         facesTmp[4] = cone[2]; facesTmp[5] = cone[0];
30*9a5b13a2SMatthew G Knepley         *faces = facesTmp;
31*9a5b13a2SMatthew G Knepley       }
32*9a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 3;
33*9a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 2;
349f074e33SMatthew G Knepley       break;
359f074e33SMatthew G Knepley     case 4:
36*9a5b13a2SMatthew G Knepley       /* Vertices follow right hand rule */
37*9a5b13a2SMatthew G Knepley       if (faces) {
38*9a5b13a2SMatthew G Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
39*9a5b13a2SMatthew G Knepley         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
40*9a5b13a2SMatthew G Knepley         facesTmp[4] = cone[2]; facesTmp[5] = cone[3];
41*9a5b13a2SMatthew G Knepley         facesTmp[6] = cone[3]; facesTmp[7] = cone[0];
42*9a5b13a2SMatthew G Knepley         *faces = facesTmp;
43*9a5b13a2SMatthew G Knepley       }
44*9a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 4;
45*9a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 2;
46*9a5b13a2SMatthew 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:
55*9a5b13a2SMatthew G Knepley       if (faces) {
56*9a5b13a2SMatthew G Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
57*9a5b13a2SMatthew G Knepley         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
58*9a5b13a2SMatthew G Knepley         facesTmp[4] = cone[2]; facesTmp[5] = cone[0];
59*9a5b13a2SMatthew G Knepley         *faces = facesTmp;
60*9a5b13a2SMatthew G Knepley       }
61*9a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 3;
62*9a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 2;
63*9a5b13a2SMatthew G Knepley       if (faces)    *faces            = facesTmp;
649f074e33SMatthew G Knepley       break;
659f074e33SMatthew G Knepley     case 4:
66*9a5b13a2SMatthew G Knepley       /* Vertices of first face follow right hand rule and normal points towards last vertex */
67*9a5b13a2SMatthew G Knepley       if (faces) {
68*9a5b13a2SMatthew G Knepley         facesTmp[0] = cone[0]; facesTmp[1]  = cone[2]; facesTmp[2]  = cone[1];
69*9a5b13a2SMatthew G Knepley         facesTmp[3] = cone[0]; facesTmp[4]  = cone[1]; facesTmp[5]  = cone[3];
70*9a5b13a2SMatthew G Knepley         facesTmp[6] = cone[0]; facesTmp[7]  = cone[3]; facesTmp[8]  = cone[2];
71*9a5b13a2SMatthew G Knepley         facesTmp[9] = cone[1]; facesTmp[10] = cone[2]; facesTmp[11] = cone[3];
72*9a5b13a2SMatthew G Knepley         *faces = facesTmp;
73*9a5b13a2SMatthew G Knepley       }
74*9a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 4;
75*9a5b13a2SMatthew G Knepley       if (faceSize) *faceSize         = 3;
76*9a5b13a2SMatthew G Knepley       if (faces)    *faces            = facesTmp;
77*9a5b13a2SMatthew G Knepley       break;
78*9a5b13a2SMatthew G Knepley     case 8:
79*9a5b13a2SMatthew G Knepley       if (faces) {
80*9a5b13a2SMatthew G Knepley         facesTmp[0]  = cone[0]; facesTmp[1]  = cone[3]; facesTmp[2]  = cone[2]; facesTmp[3]  = cone[1];
81*9a5b13a2SMatthew G Knepley         facesTmp[4]  = cone[4]; facesTmp[5]  = cone[5]; facesTmp[6]  = cone[3]; facesTmp[7]  = cone[7];
82*9a5b13a2SMatthew G Knepley         facesTmp[8]  = cone[0]; facesTmp[9]  = cone[1]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4];
83*9a5b13a2SMatthew G Knepley         facesTmp[12] = cone[2]; facesTmp[13] = cone[3]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6];
84*9a5b13a2SMatthew G Knepley         facesTmp[16] = cone[1]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5];
85*9a5b13a2SMatthew G Knepley         facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[3];
86*9a5b13a2SMatthew G Knepley         *faces = facesTmp;
87*9a5b13a2SMatthew G Knepley       }
88*9a5b13a2SMatthew G Knepley       if (numFaces) *numFaces         = 6;
89*9a5b13a2SMatthew 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   }
98*9a5b13a2SMatthew 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__
103*9a5b13a2SMatthew G Knepley #define __FUNCT__ "DMPlexInterpolateFaces_Internal"
104*9a5b13a2SMatthew G Knepley /* This interpolates faces for cells at some stratum */
105*9a5b13a2SMatthew G Knepley static PetscErrorCode DMPlexInterpolateFaces_Internal(DM dm, PetscInt cellDepth, DM idm)
1069f074e33SMatthew G Knepley {
107*9a5b13a2SMatthew G Knepley   PetscHashIJKL  faceTable;
108*9a5b13a2SMatthew G Knepley   PetscInt      *pStart, *pEnd;
109*9a5b13a2SMatthew G Knepley   PetscInt       cellDim, depth, faceDepth = cellDepth, numPoints = 0, faceSizeAll = 0, face, c, d;
1109f074e33SMatthew G Knepley   PetscErrorCode ierr;
1119f074e33SMatthew G Knepley 
1129f074e33SMatthew G Knepley   PetscFunctionBegin;
113*9a5b13a2SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &cellDim);CHKERRQ(ierr);
114*9a5b13a2SMatthew G Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
115*9a5b13a2SMatthew G Knepley   ++depth;
116*9a5b13a2SMatthew G Knepley   ++cellDepth;
117*9a5b13a2SMatthew G Knepley   cellDim -= depth - cellDepth;
118*9a5b13a2SMatthew G Knepley   ierr = PetscMalloc2(depth+1,PetscInt,&pStart,depth+1,PetscInt,&pEnd);CHKERRQ(ierr);
119*9a5b13a2SMatthew G Knepley   for (d = depth-1; d >= faceDepth; --d) {
120*9a5b13a2SMatthew G Knepley     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d+1], &pEnd[d+1]);CHKERRQ(ierr);
1219f074e33SMatthew G Knepley   }
122*9a5b13a2SMatthew G Knepley   ierr = DMPlexGetDepthStratum(dm, -1, NULL, &pStart[faceDepth]);CHKERRQ(ierr);
123*9a5b13a2SMatthew G Knepley   pEnd[faceDepth] = pStart[faceDepth];
124*9a5b13a2SMatthew G Knepley   for (d = faceDepth-1; d >= 0; --d) {
125*9a5b13a2SMatthew G Knepley     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
1269f074e33SMatthew G Knepley   }
127*9a5b13a2SMatthew G Knepley   if (pEnd[cellDepth] > pStart[cellDepth]) {ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr);}
128*9a5b13a2SMatthew 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);
129*9a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr);
130*9a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLSetMultivalued(faceTable, PETSC_FALSE);CHKERRQ(ierr);
131*9a5b13a2SMatthew G Knepley   for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) {
1329f074e33SMatthew G Knepley     const PetscInt *cellFaces;
133*9a5b13a2SMatthew G Knepley     PetscInt        numCellFaces, faceSize, cf, f;
1349f074e33SMatthew G Knepley 
135*9a5b13a2SMatthew G Knepley     ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
136*9a5b13a2SMatthew G Knepley     if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll);
1379f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
138*9a5b13a2SMatthew G Knepley       const PetscInt  *cellFace = &cellFaces[cf*faceSize];
139*9a5b13a2SMatthew G Knepley       PetscHashIJKLKey key;
1409f074e33SMatthew G Knepley 
141*9a5b13a2SMatthew G Knepley       if (faceSize == 2) {
142*9a5b13a2SMatthew G Knepley         key.i = PetscMin(cellFace[0], cellFace[1]);
143*9a5b13a2SMatthew G Knepley         key.j = PetscMax(cellFace[0], cellFace[1]);
144*9a5b13a2SMatthew G Knepley       } else {
145*9a5b13a2SMatthew G Knepley         key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0;
146*9a5b13a2SMatthew G Knepley         ierr = PetscSortInt(faceSize, (PetscInt *) &key);
1479f074e33SMatthew G Knepley       }
148*9a5b13a2SMatthew G Knepley       ierr  = PetscHashIJKLGet(faceTable, key, &f);CHKERRQ(ierr);
149*9a5b13a2SMatthew G Knepley       if (f < 0) {
150*9a5b13a2SMatthew G Knepley         ierr = PetscHashIJKLAdd(faceTable, key, face);CHKERRQ(ierr);
151*9a5b13a2SMatthew G Knepley         f    = face++;
152*9a5b13a2SMatthew G Knepley       }
153*9a5b13a2SMatthew G Knepley     }
154*9a5b13a2SMatthew G Knepley   }
155*9a5b13a2SMatthew G Knepley   pEnd[faceDepth] = face;
156*9a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr);
157*9a5b13a2SMatthew G Knepley   /* Count new points */
158*9a5b13a2SMatthew G Knepley   for (d = 0; d <= depth; ++d) {
159*9a5b13a2SMatthew G Knepley     numPoints += pEnd[d]-pStart[d];
160*9a5b13a2SMatthew G Knepley   }
161*9a5b13a2SMatthew G Knepley   ierr = DMPlexSetChart(idm, 0, numPoints);CHKERRQ(ierr);
162*9a5b13a2SMatthew G Knepley   /* Set cone sizes */
163*9a5b13a2SMatthew G Knepley   for (d = 0; d <= depth; ++d) {
164*9a5b13a2SMatthew G Knepley     PetscInt coneSize, p;
1659f074e33SMatthew G Knepley 
166*9a5b13a2SMatthew G Knepley     if (d == faceDepth) {
167*9a5b13a2SMatthew G Knepley       for (p = pStart[d]; p < pEnd[d]; ++p) {
168*9a5b13a2SMatthew G Knepley         /* I see no way to do this if we admit faces of different shapes */
169*9a5b13a2SMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr);
170*9a5b13a2SMatthew G Knepley       }
171*9a5b13a2SMatthew G Knepley     } else {
172*9a5b13a2SMatthew G Knepley       for (p = pStart[d]; p < pEnd[d]; ++p) {
173*9a5b13a2SMatthew G Knepley         ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
174*9a5b13a2SMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr);
1759f074e33SMatthew G Knepley       }
1769f074e33SMatthew G Knepley     }
1779f074e33SMatthew G Knepley   }
1789f074e33SMatthew G Knepley   ierr = DMSetUp(idm);CHKERRQ(ierr);
1799f074e33SMatthew G Knepley   /* Get face cones from subsets of cell vertices */
180*9a5b13a2SMatthew 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);
181*9a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr);
182*9a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLSetMultivalued(faceTable, PETSC_FALSE);CHKERRQ(ierr);
183*9a5b13a2SMatthew G Knepley   for (d = depth; d > cellDepth; --d) {
184*9a5b13a2SMatthew G Knepley     const PetscInt *cone;
185*9a5b13a2SMatthew G Knepley     PetscInt        p;
186*9a5b13a2SMatthew G Knepley 
187*9a5b13a2SMatthew G Knepley     for (p = pStart[d]; p < pEnd[d]; ++p) {
188*9a5b13a2SMatthew G Knepley       ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
189*9a5b13a2SMatthew G Knepley       ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr);
190*9a5b13a2SMatthew G Knepley       ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr);
191*9a5b13a2SMatthew G Knepley       ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr);
1929f074e33SMatthew G Knepley     }
193*9a5b13a2SMatthew G Knepley   }
194*9a5b13a2SMatthew G Knepley   for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) {
1959f074e33SMatthew G Knepley     const PetscInt *cellFaces;
196*9a5b13a2SMatthew G Knepley     PetscInt        numCellFaces, faceSize, cf, f;
1979f074e33SMatthew G Knepley 
198*9a5b13a2SMatthew G Knepley     ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
199*9a5b13a2SMatthew G Knepley     if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll);
2009f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
201*9a5b13a2SMatthew G Knepley       const PetscInt  *cellFace = &cellFaces[cf*faceSize];
202*9a5b13a2SMatthew G Knepley       PetscHashIJKLKey key;
2039f074e33SMatthew G Knepley 
204*9a5b13a2SMatthew G Knepley       if (faceSize == 2) {
205*9a5b13a2SMatthew G Knepley         key.i = PetscMin(cellFace[0], cellFace[1]);
206*9a5b13a2SMatthew G Knepley         key.j = PetscMax(cellFace[0], cellFace[1]);
207*9a5b13a2SMatthew G Knepley       } else {
208*9a5b13a2SMatthew G Knepley         key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0;
209*9a5b13a2SMatthew G Knepley         ierr = PetscSortInt(faceSize, (PetscInt *) &key);
2109f074e33SMatthew G Knepley       }
211*9a5b13a2SMatthew G Knepley       ierr  = PetscHashIJKLGet(faceTable, key, &f);CHKERRQ(ierr);
212*9a5b13a2SMatthew G Knepley       if (f < 0) {
213*9a5b13a2SMatthew G Knepley         ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr);
214*9a5b13a2SMatthew G Knepley         ierr = PetscHashIJKLAdd(faceTable, key, face);CHKERRQ(ierr);
215*9a5b13a2SMatthew G Knepley         f    = face++;
2169f074e33SMatthew G Knepley         ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr);
217*9a5b13a2SMatthew G Knepley       } else {
218*9a5b13a2SMatthew G Knepley         const PetscInt *cone;
219*9a5b13a2SMatthew G Knepley         PetscInt        coneSize, ornt, i, j;
2209f074e33SMatthew G Knepley 
221*9a5b13a2SMatthew G Knepley         ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr);
222*9a5b13a2SMatthew G Knepley         /* Orient face */
2239f074e33SMatthew G Knepley         ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr);
2249f074e33SMatthew G Knepley         ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr);
225*9a5b13a2SMatthew 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);
226*9a5b13a2SMatthew G Knepley         /* - First find the initial vertex */
227*9a5b13a2SMatthew G Knepley         for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break;
228*9a5b13a2SMatthew G Knepley         /* - Try forward comparison */
229*9a5b13a2SMatthew G Knepley         for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break;
230*9a5b13a2SMatthew G Knepley         if (j == faceSize) {
231*9a5b13a2SMatthew G Knepley           if ((faceSize == 2) && (i == 1)) ornt = -2;
232*9a5b13a2SMatthew G Knepley           else                             ornt = i;
233*9a5b13a2SMatthew G Knepley         } else {
234*9a5b13a2SMatthew G Knepley           /* - Try backward comparison */
235*9a5b13a2SMatthew G Knepley           for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break;
236*9a5b13a2SMatthew G Knepley           if (j == faceSize) ornt = -(i+1);
237*9a5b13a2SMatthew G Knepley           else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation");
238*9a5b13a2SMatthew G Knepley         }
239*9a5b13a2SMatthew G Knepley         ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr);
2409f074e33SMatthew G Knepley       }
2419f074e33SMatthew G Knepley     }
2429f074e33SMatthew G Knepley   }
243*9a5b13a2SMatthew 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]);
244*9a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr);
245*9a5b13a2SMatthew G Knepley   ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr);
246*9a5b13a2SMatthew G Knepley   ierr = DMPlexStratify(idm);CHKERRQ(ierr);
2479f074e33SMatthew G Knepley   PetscFunctionReturn(0);
2489f074e33SMatthew G Knepley }
2499f074e33SMatthew G Knepley 
2509f074e33SMatthew G Knepley #undef __FUNCT__
2519f074e33SMatthew G Knepley #define __FUNCT__ "DMPlexInterpolate"
2529f074e33SMatthew G Knepley PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt)
2539f074e33SMatthew G Knepley {
254*9a5b13a2SMatthew G Knepley   DM             idm, odm = dm;
255*9a5b13a2SMatthew G Knepley   PetscInt       dim, d;
2569f074e33SMatthew G Knepley   PetscErrorCode ierr;
2579f074e33SMatthew G Knepley 
2589f074e33SMatthew G Knepley   PetscFunctionBegin;
2599f074e33SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
260*9a5b13a2SMatthew G Knepley   for (d = 1; d < dim; ++d) {
261*9a5b13a2SMatthew G Knepley     /* Create interpolated mesh */
262*9a5b13a2SMatthew G Knepley     ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr);
263*9a5b13a2SMatthew G Knepley     ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr);
264*9a5b13a2SMatthew G Knepley     ierr = DMPlexSetDimension(idm, dim);CHKERRQ(ierr);
265*9a5b13a2SMatthew G Knepley     ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr);
266*9a5b13a2SMatthew G Knepley     if (dim == 3) {ierr = DMView(idm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);}
267*9a5b13a2SMatthew G Knepley     if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);}
268*9a5b13a2SMatthew G Knepley     odm  = idm;
2699f074e33SMatthew G Knepley   }
270*9a5b13a2SMatthew G Knepley   *dmInt = idm;
2719f074e33SMatthew G Knepley   PetscFunctionReturn(0);
2729f074e33SMatthew G Knepley }
27307b0a1fcSMatthew G Knepley 
27407b0a1fcSMatthew G Knepley #undef __FUNCT__
27507b0a1fcSMatthew G Knepley #define __FUNCT__ "DMPlexCopyCoordinates"
27607b0a1fcSMatthew G Knepley PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB)
27707b0a1fcSMatthew G Knepley {
27807b0a1fcSMatthew G Knepley   Vec            coordinatesA, coordinatesB;
27907b0a1fcSMatthew G Knepley   PetscSection   coordSectionA, coordSectionB;
28007b0a1fcSMatthew G Knepley   PetscScalar   *coordsA, *coordsB;
28107b0a1fcSMatthew G Knepley   PetscInt       spaceDim, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d;
28207b0a1fcSMatthew G Knepley   PetscErrorCode ierr;
28307b0a1fcSMatthew G Knepley 
28407b0a1fcSMatthew G Knepley   PetscFunctionBegin;
28507b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr);
28607b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr);
28707b0a1fcSMatthew 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);
28807b0a1fcSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr);
28907b0a1fcSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr);
29007b0a1fcSMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr);
29107b0a1fcSMatthew G Knepley   ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr);
29207b0a1fcSMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr);
29307b0a1fcSMatthew G Knepley   ierr = PetscSectionSetChart(coordSectionB, vStartB, vEndB);CHKERRQ(ierr);
29407b0a1fcSMatthew G Knepley   for (v = vStartB; v < vEndB; ++v) {
29507b0a1fcSMatthew G Knepley     ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr);
29607b0a1fcSMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr);
29707b0a1fcSMatthew G Knepley   }
29807b0a1fcSMatthew G Knepley   ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr);
29907b0a1fcSMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr);
30007b0a1fcSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr);
30107b0a1fcSMatthew G Knepley   ierr = VecCreate(PetscObjectComm((PetscObject) dmB), &coordinatesB);CHKERRQ(ierr);
30207b0a1fcSMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr);
30307b0a1fcSMatthew G Knepley   ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr);
30407b0a1fcSMatthew G Knepley   ierr = VecSetFromOptions(coordinatesB);CHKERRQ(ierr);
30507b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr);
30607b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr);
30707b0a1fcSMatthew G Knepley   for (v = 0; v < vEndB-vStartB; ++v) {
30807b0a1fcSMatthew G Knepley     for (d = 0; d < spaceDim; ++d) {
30907b0a1fcSMatthew G Knepley       coordsB[v*spaceDim+d] = coordsA[v*spaceDim+d];
31007b0a1fcSMatthew G Knepley     }
31107b0a1fcSMatthew G Knepley   }
31207b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr);
31307b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr);
31407b0a1fcSMatthew G Knepley   ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr);
31507b0a1fcSMatthew G Knepley   ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr);
31607b0a1fcSMatthew G Knepley   PetscFunctionReturn(0);
31707b0a1fcSMatthew G Knepley }
318