xref: /petsc/src/dm/impls/plex/plexinterpolate.c (revision 9f074e33fb8b4d6fdad352a0506a3b1288f5220f)
1*9f074e33SMatthew G Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2*9f074e33SMatthew G Knepley #include <../src/sys/utils/hash.h>
3*9f074e33SMatthew G Knepley 
4*9f074e33SMatthew G Knepley #undef __FUNCT__
5*9f074e33SMatthew G Knepley #define __FUNCT__ "DMPlexGetFaces"
6*9f074e33SMatthew G Knepley /*
7*9f074e33SMatthew G Knepley   DMPlexGetFaces -
8*9f074e33SMatthew G Knepley 
9*9f074e33SMatthew G Knepley   Note: This will only work for cell-vertex meshes.
10*9f074e33SMatthew G Knepley */
11*9f074e33SMatthew G Knepley static PetscErrorCode DMPlexGetFaces(DM dm, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[])
12*9f074e33SMatthew G Knepley {
13*9f074e33SMatthew G Knepley   DM_Plex        *mesh = (DM_Plex*) dm->data;
14*9f074e33SMatthew G Knepley   const PetscInt *cone = NULL;
15*9f074e33SMatthew G Knepley   PetscInt        depth = 0, dim, coneSize;
16*9f074e33SMatthew G Knepley   PetscErrorCode  ierr;
17*9f074e33SMatthew G Knepley 
18*9f074e33SMatthew G Knepley   PetscFunctionBegin;
19*9f074e33SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
20*9f074e33SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
21*9f074e33SMatthew G Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
22*9f074e33SMatthew G Knepley   if (depth > 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Faces can only be returned for cell-vertex meshes.");
23*9f074e33SMatthew G Knepley   if (!mesh->facesTmp) {ierr = PetscMalloc(PetscSqr(PetscMax(mesh->maxConeSize, mesh->maxSupportSize)) * sizeof(PetscInt), &mesh->facesTmp);CHKERRQ(ierr);}
24*9f074e33SMatthew G Knepley   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
25*9f074e33SMatthew G Knepley   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
26*9f074e33SMatthew G Knepley   switch (dim) {
27*9f074e33SMatthew G Knepley   case 2:
28*9f074e33SMatthew G Knepley     switch (coneSize) {
29*9f074e33SMatthew G Knepley     case 3:
30*9f074e33SMatthew G Knepley       mesh->facesTmp[0] = cone[0]; mesh->facesTmp[1] = cone[1];
31*9f074e33SMatthew G Knepley       mesh->facesTmp[2] = cone[1]; mesh->facesTmp[3] = cone[2];
32*9f074e33SMatthew G Knepley       mesh->facesTmp[4] = cone[2]; mesh->facesTmp[5] = cone[0];
33*9f074e33SMatthew G Knepley       *numFaces         = 3;
34*9f074e33SMatthew G Knepley       *faceSize         = 2;
35*9f074e33SMatthew G Knepley       *faces            = mesh->facesTmp;
36*9f074e33SMatthew G Knepley       break;
37*9f074e33SMatthew G Knepley     case 4:
38*9f074e33SMatthew G Knepley       mesh->facesTmp[0] = cone[0]; mesh->facesTmp[1] = cone[1];
39*9f074e33SMatthew G Knepley       mesh->facesTmp[2] = cone[1]; mesh->facesTmp[3] = cone[2];
40*9f074e33SMatthew G Knepley       mesh->facesTmp[4] = cone[2]; mesh->facesTmp[5] = cone[3];
41*9f074e33SMatthew G Knepley       mesh->facesTmp[6] = cone[3]; mesh->facesTmp[7] = cone[0];
42*9f074e33SMatthew G Knepley       *numFaces         = 4;
43*9f074e33SMatthew G Knepley       *faceSize         = 2;
44*9f074e33SMatthew G Knepley       *faces            = mesh->facesTmp;
45*9f074e33SMatthew G Knepley       break;
46*9f074e33SMatthew G Knepley     default:
47*9f074e33SMatthew G Knepley       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
48*9f074e33SMatthew G Knepley     }
49*9f074e33SMatthew G Knepley     break;
50*9f074e33SMatthew G Knepley   case 3:
51*9f074e33SMatthew G Knepley     switch (coneSize) {
52*9f074e33SMatthew G Knepley     case 3:
53*9f074e33SMatthew G Knepley       mesh->facesTmp[0] = cone[0]; mesh->facesTmp[1] = cone[1];
54*9f074e33SMatthew G Knepley       mesh->facesTmp[2] = cone[1]; mesh->facesTmp[3] = cone[2];
55*9f074e33SMatthew G Knepley       mesh->facesTmp[4] = cone[2]; mesh->facesTmp[5] = cone[0];
56*9f074e33SMatthew G Knepley       *numFaces         = 3;
57*9f074e33SMatthew G Knepley       *faceSize         = 2;
58*9f074e33SMatthew G Knepley       *faces            = mesh->facesTmp;
59*9f074e33SMatthew G Knepley       break;
60*9f074e33SMatthew G Knepley     case 4:
61*9f074e33SMatthew G Knepley       mesh->facesTmp[0] = cone[0]; mesh->facesTmp[1]  = cone[1]; mesh->facesTmp[2]  = cone[2];
62*9f074e33SMatthew G Knepley       mesh->facesTmp[3] = cone[0]; mesh->facesTmp[4]  = cone[2]; mesh->facesTmp[5]  = cone[3];
63*9f074e33SMatthew G Knepley       mesh->facesTmp[6] = cone[0]; mesh->facesTmp[7]  = cone[3]; mesh->facesTmp[8]  = cone[1];
64*9f074e33SMatthew G Knepley       mesh->facesTmp[9] = cone[1]; mesh->facesTmp[10] = cone[3]; mesh->facesTmp[11] = cone[2];
65*9f074e33SMatthew G Knepley       *numFaces         = 4;
66*9f074e33SMatthew G Knepley       *faceSize         = 3;
67*9f074e33SMatthew G Knepley       *faces            = mesh->facesTmp;
68*9f074e33SMatthew G Knepley       break;
69*9f074e33SMatthew G Knepley     default:
70*9f074e33SMatthew G Knepley       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
71*9f074e33SMatthew G Knepley     }
72*9f074e33SMatthew G Knepley     break;
73*9f074e33SMatthew G Knepley   default:
74*9f074e33SMatthew G Knepley     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim);
75*9f074e33SMatthew G Knepley   }
76*9f074e33SMatthew G Knepley   PetscFunctionReturn(0);
77*9f074e33SMatthew G Knepley }
78*9f074e33SMatthew G Knepley 
79*9f074e33SMatthew G Knepley #undef __FUNCT__
80*9f074e33SMatthew G Knepley #define __FUNCT__ "DMPlexInterpolate_2D"
81*9f074e33SMatthew G Knepley static PetscErrorCode DMPlexInterpolate_2D(DM dm, DM *dmInt)
82*9f074e33SMatthew G Knepley {
83*9f074e33SMatthew G Knepley   DM             idm;
84*9f074e33SMatthew G Knepley   DM_Plex       *mesh;
85*9f074e33SMatthew G Knepley   PetscHashIJ    edgeTable;
86*9f074e33SMatthew G Knepley   PetscInt      *off;
87*9f074e33SMatthew G Knepley   PetscInt       dim, numCells, cStart, cEnd, c, numVertices, vStart, vEnd;
88*9f074e33SMatthew G Knepley   PetscInt       numEdges, firstEdge, edge, e;
89*9f074e33SMatthew G Knepley   PetscErrorCode ierr;
90*9f074e33SMatthew G Knepley 
91*9f074e33SMatthew G Knepley   PetscFunctionBegin;
92*9f074e33SMatthew G Knepley   ierr        = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
93*9f074e33SMatthew G Knepley   ierr        = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
94*9f074e33SMatthew G Knepley   ierr        = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
95*9f074e33SMatthew G Knepley   numCells    = cEnd - cStart;
96*9f074e33SMatthew G Knepley   numVertices = vEnd - vStart;
97*9f074e33SMatthew G Knepley   firstEdge   = numCells + numVertices;
98*9f074e33SMatthew G Knepley   numEdges    = 0;
99*9f074e33SMatthew G Knepley   /* Count edges using algorithm from CreateNeighborCSR */
100*9f074e33SMatthew G Knepley   ierr = DMPlexCreateNeighborCSR(dm, NULL, &off, NULL);CHKERRQ(ierr);
101*9f074e33SMatthew G Knepley   if (off) {
102*9f074e33SMatthew G Knepley     PetscInt numCorners = 0;
103*9f074e33SMatthew G Knepley 
104*9f074e33SMatthew G Knepley     numEdges = off[numCells]/2;
105*9f074e33SMatthew G Knepley #if 0
106*9f074e33SMatthew G Knepley     /* Account for boundary edges: \sum_c 3 - neighbors = 3*numCells - totalNeighbors */
107*9f074e33SMatthew G Knepley     numEdges += 3*numCells - off[numCells];
108*9f074e33SMatthew G Knepley #else
109*9f074e33SMatthew G Knepley     /* Account for boundary edges: \sum_c #faces - #neighbors = \sum_c #cellVertices - #neighbors = totalCorners - totalNeighbors */
110*9f074e33SMatthew G Knepley     for (c = cStart; c < cEnd; ++c) {
111*9f074e33SMatthew G Knepley       PetscInt coneSize;
112*9f074e33SMatthew G Knepley 
113*9f074e33SMatthew G Knepley       ierr        = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
114*9f074e33SMatthew G Knepley       numCorners += coneSize;
115*9f074e33SMatthew G Knepley     }
116*9f074e33SMatthew G Knepley     numEdges += numCorners - off[numCells];
117*9f074e33SMatthew G Knepley #endif
118*9f074e33SMatthew G Knepley   }
119*9f074e33SMatthew G Knepley #if 0
120*9f074e33SMatthew G Knepley   /* Check Euler characteristic V - E + F = 1 */
121*9f074e33SMatthew G Knepley   if (numVertices && (numVertices-numEdges+numCells != 1)) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Euler characteristic of mesh is %d  != 1", numVertices-numEdges+numCells);
122*9f074e33SMatthew G Knepley #endif
123*9f074e33SMatthew G Knepley   /* Create interpolated mesh */
124*9f074e33SMatthew G Knepley   ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr);
125*9f074e33SMatthew G Knepley   ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr);
126*9f074e33SMatthew G Knepley   ierr = DMPlexSetDimension(idm, dim);CHKERRQ(ierr);
127*9f074e33SMatthew G Knepley   ierr = DMPlexSetChart(idm, 0, numCells+numVertices+numEdges);CHKERRQ(ierr);
128*9f074e33SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
129*9f074e33SMatthew G Knepley     PetscInt numCorners;
130*9f074e33SMatthew G Knepley 
131*9f074e33SMatthew G Knepley     ierr = DMPlexGetConeSize(dm, c, &numCorners);CHKERRQ(ierr);
132*9f074e33SMatthew G Knepley     ierr = DMPlexSetConeSize(idm, c, numCorners);CHKERRQ(ierr);
133*9f074e33SMatthew G Knepley   }
134*9f074e33SMatthew G Knepley   for (e = firstEdge; e < firstEdge+numEdges; ++e) {
135*9f074e33SMatthew G Knepley     ierr = DMPlexSetConeSize(idm, e, 2);CHKERRQ(ierr);
136*9f074e33SMatthew G Knepley   }
137*9f074e33SMatthew G Knepley   ierr = DMSetUp(idm);CHKERRQ(ierr);
138*9f074e33SMatthew G Knepley   /* Get edge cones from subsets of cell vertices */
139*9f074e33SMatthew G Knepley   ierr = PetscHashIJCreate(&edgeTable);CHKERRQ(ierr);
140*9f074e33SMatthew G Knepley   ierr = PetscHashIJSetMultivalued(edgeTable, PETSC_FALSE);CHKERRQ(ierr);
141*9f074e33SMatthew G Knepley 
142*9f074e33SMatthew G Knepley   for (c = 0, edge = firstEdge; c < numCells; ++c) {
143*9f074e33SMatthew G Knepley     const PetscInt *cellFaces;
144*9f074e33SMatthew G Knepley     PetscInt        numCellFaces, faceSize, cf;
145*9f074e33SMatthew G Knepley 
146*9f074e33SMatthew G Knepley     ierr = DMPlexGetFaces(dm, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
147*9f074e33SMatthew G Knepley     if (faceSize != 2) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Triangles cannot have face of size %D", faceSize);
148*9f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
149*9f074e33SMatthew G Knepley #if 1
150*9f074e33SMatthew G Knepley       PetscHashIJKey key;
151*9f074e33SMatthew G Knepley 
152*9f074e33SMatthew G Knepley       key.i = PetscMin(cellFaces[cf*faceSize+0], cellFaces[cf*faceSize+1]);
153*9f074e33SMatthew G Knepley       key.j = PetscMax(cellFaces[cf*faceSize+0], cellFaces[cf*faceSize+1]);
154*9f074e33SMatthew G Knepley       ierr  = PetscHashIJGet(edgeTable, key, &e);CHKERRQ(ierr);
155*9f074e33SMatthew G Knepley       if (e < 0) {
156*9f074e33SMatthew G Knepley         ierr = DMPlexSetCone(idm, edge, &cellFaces[cf*faceSize]);CHKERRQ(ierr);
157*9f074e33SMatthew G Knepley         ierr = PetscHashIJAdd(edgeTable, key, edge);CHKERRQ(ierr);
158*9f074e33SMatthew G Knepley         e    = edge++;
159*9f074e33SMatthew G Knepley       }
160*9f074e33SMatthew G Knepley #else
161*9f074e33SMatthew G Knepley       PetscBool found = PETSC_FALSE;
162*9f074e33SMatthew G Knepley 
163*9f074e33SMatthew G Knepley       /* TODO Need join of vertices to check for existence of edges, which needs support (could set edge support), so just brute force for now */
164*9f074e33SMatthew G Knepley       for (e = firstEdge; e < edge; ++e) {
165*9f074e33SMatthew G Knepley         const PetscInt *cone;
166*9f074e33SMatthew G Knepley 
167*9f074e33SMatthew G Knepley         ierr = DMPlexGetCone(idm, e, &cone);CHKERRQ(ierr);
168*9f074e33SMatthew G Knepley         if (((cellFaces[cf*faceSize+0] == cone[0]) && (cellFaces[cf*faceSize+1] == cone[1])) ||
169*9f074e33SMatthew G Knepley             ((cellFaces[cf*faceSize+0] == cone[1]) && (cellFaces[cf*faceSize+1] == cone[0]))) {
170*9f074e33SMatthew G Knepley           found = PETSC_TRUE;
171*9f074e33SMatthew G Knepley           break;
172*9f074e33SMatthew G Knepley         }
173*9f074e33SMatthew G Knepley       }
174*9f074e33SMatthew G Knepley       if (!found) {
175*9f074e33SMatthew G Knepley         ierr = DMPlexSetCone(idm, edge, &cellFaces[cf*faceSize]);CHKERRQ(ierr);
176*9f074e33SMatthew G Knepley         ++edge;
177*9f074e33SMatthew G Knepley       }
178*9f074e33SMatthew G Knepley #endif
179*9f074e33SMatthew G Knepley       ierr = DMPlexInsertCone(idm, c, cf, e);CHKERRQ(ierr);
180*9f074e33SMatthew G Knepley     }
181*9f074e33SMatthew G Knepley   }
182*9f074e33SMatthew G Knepley   if (edge != firstEdge+numEdges) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid number of edges %D should be %D", edge-firstEdge, numEdges);
183*9f074e33SMatthew G Knepley   ierr = PetscHashIJDestroy(&edgeTable);CHKERRQ(ierr);
184*9f074e33SMatthew G Knepley   ierr = PetscFree(off);CHKERRQ(ierr);
185*9f074e33SMatthew G Knepley   ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr);
186*9f074e33SMatthew G Knepley   ierr = DMPlexStratify(idm);CHKERRQ(ierr);
187*9f074e33SMatthew G Knepley   mesh = (DM_Plex*) (idm)->data;
188*9f074e33SMatthew G Knepley   /* Orient edges */
189*9f074e33SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
190*9f074e33SMatthew G Knepley     const PetscInt *cone = NULL, *cellFaces;
191*9f074e33SMatthew G Knepley     PetscInt        coneSize, coff, numCellFaces, faceSize, cf;
192*9f074e33SMatthew G Knepley 
193*9f074e33SMatthew G Knepley     ierr = DMPlexGetConeSize(idm, c, &coneSize);CHKERRQ(ierr);
194*9f074e33SMatthew G Knepley     ierr = DMPlexGetCone(idm, c, &cone);CHKERRQ(ierr);
195*9f074e33SMatthew G Knepley     ierr = PetscSectionGetOffset(mesh->coneSection, c, &coff);CHKERRQ(ierr);
196*9f074e33SMatthew G Knepley     ierr = DMPlexGetFaces(dm, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
197*9f074e33SMatthew G Knepley     if (coneSize != numCellFaces) SETERRQ3(PetscObjectComm((PetscObject)idm), PETSC_ERR_PLIB, "Invalid number of edges %D for cell %D should be %D", coneSize, c, numCellFaces);
198*9f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
199*9f074e33SMatthew G Knepley       const PetscInt *econe = NULL;
200*9f074e33SMatthew G Knepley       PetscInt        esize;
201*9f074e33SMatthew G Knepley 
202*9f074e33SMatthew G Knepley       ierr = DMPlexGetConeSize(idm, cone[cf], &esize);CHKERRQ(ierr);
203*9f074e33SMatthew G Knepley       ierr = DMPlexGetCone(idm, cone[cf], &econe);CHKERRQ(ierr);
204*9f074e33SMatthew G Knepley       if (esize != 2) SETERRQ2(PetscObjectComm((PetscObject)idm), PETSC_ERR_PLIB, "Invalid number of edge endpoints %D for edge %D should be 2", esize, cone[cf]);
205*9f074e33SMatthew G Knepley       if ((cellFaces[cf*faceSize+0] == econe[0]) && (cellFaces[cf*faceSize+1] == econe[1])) {
206*9f074e33SMatthew G Knepley         /* Correctly oriented */
207*9f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = 0;
208*9f074e33SMatthew G Knepley       } else if ((cellFaces[cf*faceSize+0] == econe[1]) && (cellFaces[cf*faceSize+1] == econe[0])) {
209*9f074e33SMatthew G Knepley         /* Start at index 1, and reverse orientation */
210*9f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = -(1+1);
211*9f074e33SMatthew G Knepley       }
212*9f074e33SMatthew G Knepley     }
213*9f074e33SMatthew G Knepley   }
214*9f074e33SMatthew G Knepley   *dmInt = idm;
215*9f074e33SMatthew G Knepley   PetscFunctionReturn(0);
216*9f074e33SMatthew G Knepley }
217*9f074e33SMatthew G Knepley 
218*9f074e33SMatthew G Knepley #undef __FUNCT__
219*9f074e33SMatthew G Knepley #define __FUNCT__ "DMPlexInterpolate_3D"
220*9f074e33SMatthew G Knepley static PetscErrorCode DMPlexInterpolate_3D(DM dm, DM *dmInt)
221*9f074e33SMatthew G Knepley {
222*9f074e33SMatthew G Knepley   DM             idm, fdm;
223*9f074e33SMatthew G Knepley   DM_Plex       *mesh;
224*9f074e33SMatthew G Knepley   PetscInt      *off;
225*9f074e33SMatthew G Knepley   const PetscInt numCorners = 4;
226*9f074e33SMatthew G Knepley   PetscInt       dim, numCells, cStart, cEnd, c, numVertices, vStart, vEnd;
227*9f074e33SMatthew G Knepley   PetscInt       numFaces, firstFace, face, f, numEdges, firstEdge, edge, e;
228*9f074e33SMatthew G Knepley   PetscErrorCode ierr;
229*9f074e33SMatthew G Knepley 
230*9f074e33SMatthew G Knepley   PetscFunctionBegin;
231*9f074e33SMatthew G Knepley   {
232*9f074e33SMatthew G Knepley     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr);
233*9f074e33SMatthew G Knepley     ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
234*9f074e33SMatthew G Knepley     ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
235*9f074e33SMatthew G Knepley   }
236*9f074e33SMatthew G Knepley   ierr        = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
237*9f074e33SMatthew G Knepley   ierr        = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
238*9f074e33SMatthew G Knepley   ierr        = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
239*9f074e33SMatthew G Knepley   numCells    = cEnd - cStart;
240*9f074e33SMatthew G Knepley   numVertices = vEnd - vStart;
241*9f074e33SMatthew G Knepley   firstFace   = numCells + numVertices;
242*9f074e33SMatthew G Knepley   numFaces    = 0;
243*9f074e33SMatthew G Knepley   /* Count faces using algorithm from CreateNeighborCSR */
244*9f074e33SMatthew G Knepley   ierr = DMPlexCreateNeighborCSR(dm, NULL, &off, NULL);CHKERRQ(ierr);
245*9f074e33SMatthew G Knepley   if (off) {
246*9f074e33SMatthew G Knepley     numFaces = off[numCells]/2;
247*9f074e33SMatthew G Knepley     /* Account for boundary faces: \sum_c 4 - neighbors = 4*numCells - totalNeighbors */
248*9f074e33SMatthew G Knepley     numFaces += 4*numCells - off[numCells];
249*9f074e33SMatthew G Knepley   }
250*9f074e33SMatthew G Knepley   /* Use Euler characteristic to get edges V - E + F - C = 1 */
251*9f074e33SMatthew G Knepley   firstEdge = firstFace + numFaces;
252*9f074e33SMatthew G Knepley   numEdges  = numVertices + numFaces - numCells - 1;
253*9f074e33SMatthew G Knepley   /* Create interpolated mesh */
254*9f074e33SMatthew G Knepley   ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr);
255*9f074e33SMatthew G Knepley   ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr);
256*9f074e33SMatthew G Knepley   ierr = DMPlexSetDimension(idm, dim);CHKERRQ(ierr);
257*9f074e33SMatthew G Knepley   ierr = DMPlexSetChart(idm, 0, numCells+numVertices+numFaces+numEdges);CHKERRQ(ierr);
258*9f074e33SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
259*9f074e33SMatthew G Knepley     ierr = DMPlexSetConeSize(idm, c, numCorners);CHKERRQ(ierr);
260*9f074e33SMatthew G Knepley   }
261*9f074e33SMatthew G Knepley   for (f = firstFace; f < firstFace+numFaces; ++f) {
262*9f074e33SMatthew G Knepley     ierr = DMPlexSetConeSize(idm, f, 3);CHKERRQ(ierr);
263*9f074e33SMatthew G Knepley   }
264*9f074e33SMatthew G Knepley   for (e = firstEdge; e < firstEdge+numEdges; ++e) {
265*9f074e33SMatthew G Knepley     ierr = DMPlexSetConeSize(idm, e, 2);CHKERRQ(ierr);
266*9f074e33SMatthew G Knepley   }
267*9f074e33SMatthew G Knepley   ierr = DMSetUp(idm);CHKERRQ(ierr);
268*9f074e33SMatthew G Knepley   /* Get face cones from subsets of cell vertices */
269*9f074e33SMatthew G Knepley   ierr = DMCreate(PetscObjectComm((PetscObject)dm), &fdm);CHKERRQ(ierr);
270*9f074e33SMatthew G Knepley   ierr = DMSetType(fdm, DMPLEX);CHKERRQ(ierr);
271*9f074e33SMatthew G Knepley   ierr = DMPlexSetDimension(fdm, dim);CHKERRQ(ierr);
272*9f074e33SMatthew G Knepley   ierr = DMPlexSetChart(fdm, numCells, firstFace+numFaces);CHKERRQ(ierr);
273*9f074e33SMatthew G Knepley   for (f = firstFace; f < firstFace+numFaces; ++f) {
274*9f074e33SMatthew G Knepley     ierr = DMPlexSetConeSize(fdm, f, 3);CHKERRQ(ierr);
275*9f074e33SMatthew G Knepley   }
276*9f074e33SMatthew G Knepley   ierr = DMSetUp(fdm);CHKERRQ(ierr);
277*9f074e33SMatthew G Knepley   for (c = 0, face = firstFace; c < numCells; ++c) {
278*9f074e33SMatthew G Knepley     const PetscInt *cellFaces;
279*9f074e33SMatthew G Knepley     PetscInt        numCellFaces, faceSize, cf;
280*9f074e33SMatthew G Knepley 
281*9f074e33SMatthew G Knepley     ierr = DMPlexGetFaces(dm, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
282*9f074e33SMatthew G Knepley     if (faceSize != 3) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Tetrahedra cannot have face of size %D", faceSize);
283*9f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
284*9f074e33SMatthew G Knepley       PetscBool found = PETSC_FALSE;
285*9f074e33SMatthew G Knepley 
286*9f074e33SMatthew G Knepley       /* TODO Need join of vertices to check for existence of edges, which needs support (could set edge support), so just brute force for now */
287*9f074e33SMatthew G Knepley       for (f = firstFace; f < face; ++f) {
288*9f074e33SMatthew G Knepley         const PetscInt *cone = NULL;
289*9f074e33SMatthew G Knepley 
290*9f074e33SMatthew G Knepley         ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr);
291*9f074e33SMatthew G Knepley         if (((cellFaces[cf*faceSize+0] == cone[0]) && (cellFaces[cf*faceSize+1] == cone[1]) && (cellFaces[cf*faceSize+2] == cone[2])) ||
292*9f074e33SMatthew G Knepley             ((cellFaces[cf*faceSize+0] == cone[1]) && (cellFaces[cf*faceSize+1] == cone[2]) && (cellFaces[cf*faceSize+2] == cone[0])) ||
293*9f074e33SMatthew G Knepley             ((cellFaces[cf*faceSize+0] == cone[2]) && (cellFaces[cf*faceSize+1] == cone[0]) && (cellFaces[cf*faceSize+2] == cone[1])) ||
294*9f074e33SMatthew G Knepley             ((cellFaces[cf*faceSize+0] == cone[0]) && (cellFaces[cf*faceSize+1] == cone[2]) && (cellFaces[cf*faceSize+2] == cone[1])) ||
295*9f074e33SMatthew G Knepley             ((cellFaces[cf*faceSize+0] == cone[2]) && (cellFaces[cf*faceSize+1] == cone[1]) && (cellFaces[cf*faceSize+2] == cone[0])) ||
296*9f074e33SMatthew G Knepley             ((cellFaces[cf*faceSize+0] == cone[1]) && (cellFaces[cf*faceSize+1] == cone[0]) && (cellFaces[cf*faceSize+2] == cone[2]))) {
297*9f074e33SMatthew G Knepley           found = PETSC_TRUE;
298*9f074e33SMatthew G Knepley           break;
299*9f074e33SMatthew G Knepley         }
300*9f074e33SMatthew G Knepley       }
301*9f074e33SMatthew G Knepley       if (!found) {
302*9f074e33SMatthew G Knepley         ierr = DMPlexSetCone(idm, face, &cellFaces[cf*faceSize]);CHKERRQ(ierr);
303*9f074e33SMatthew G Knepley         /* Save the vertices for orientation calculation */
304*9f074e33SMatthew G Knepley         ierr = DMPlexSetCone(fdm, face, &cellFaces[cf*faceSize]);CHKERRQ(ierr);
305*9f074e33SMatthew G Knepley         ++face;
306*9f074e33SMatthew G Knepley       }
307*9f074e33SMatthew G Knepley       ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr);
308*9f074e33SMatthew G Knepley     }
309*9f074e33SMatthew G Knepley   }
310*9f074e33SMatthew G Knepley   if (face != firstFace+numFaces) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid number of faces %D should be %D", face-firstFace, numFaces);
311*9f074e33SMatthew G Knepley   /* Get edge cones from subsets of face vertices */
312*9f074e33SMatthew G Knepley   for (f = firstFace, edge = firstEdge; f < firstFace+numFaces; ++f) {
313*9f074e33SMatthew G Knepley     const PetscInt *cellFaces;
314*9f074e33SMatthew G Knepley     PetscInt        numCellFaces, faceSize, cf;
315*9f074e33SMatthew G Knepley 
316*9f074e33SMatthew G Knepley     ierr = DMPlexGetFaces(idm, f, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
317*9f074e33SMatthew G Knepley     if (faceSize != 2) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Triangles cannot have face of size %D", faceSize);
318*9f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
319*9f074e33SMatthew G Knepley       PetscBool found = PETSC_FALSE;
320*9f074e33SMatthew G Knepley 
321*9f074e33SMatthew G Knepley       /* TODO Need join of vertices to check for existence of edges, which needs support (could set edge support), so just brute force for now */
322*9f074e33SMatthew G Knepley       for (e = firstEdge; e < edge; ++e) {
323*9f074e33SMatthew G Knepley         const PetscInt *cone = NULL;
324*9f074e33SMatthew G Knepley 
325*9f074e33SMatthew G Knepley         ierr = DMPlexGetCone(idm, e, &cone);CHKERRQ(ierr);
326*9f074e33SMatthew G Knepley         if (((cellFaces[cf*faceSize+0] == cone[0]) && (cellFaces[cf*faceSize+1] == cone[1])) ||
327*9f074e33SMatthew G Knepley             ((cellFaces[cf*faceSize+0] == cone[1]) && (cellFaces[cf*faceSize+1] == cone[0]))) {
328*9f074e33SMatthew G Knepley           found = PETSC_TRUE;
329*9f074e33SMatthew G Knepley           break;
330*9f074e33SMatthew G Knepley         }
331*9f074e33SMatthew G Knepley       }
332*9f074e33SMatthew G Knepley       if (!found) {
333*9f074e33SMatthew G Knepley         ierr = DMPlexSetCone(idm, edge, &cellFaces[cf*faceSize]);CHKERRQ(ierr);
334*9f074e33SMatthew G Knepley         ++edge;
335*9f074e33SMatthew G Knepley       }
336*9f074e33SMatthew G Knepley       ierr = DMPlexInsertCone(idm, f, cf, e);CHKERRQ(ierr);
337*9f074e33SMatthew G Knepley     }
338*9f074e33SMatthew G Knepley   }
339*9f074e33SMatthew G Knepley   if (edge != firstEdge+numEdges) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid number of edges %D should be %D", edge-firstEdge, numEdges);
340*9f074e33SMatthew G Knepley   ierr = PetscFree(off);CHKERRQ(ierr);
341*9f074e33SMatthew G Knepley   ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr);
342*9f074e33SMatthew G Knepley   ierr = DMPlexStratify(idm);CHKERRQ(ierr);
343*9f074e33SMatthew G Knepley   mesh = (DM_Plex*) (idm)->data;
344*9f074e33SMatthew G Knepley   /* Orient edges */
345*9f074e33SMatthew G Knepley   for (f = firstFace; f < firstFace+numFaces; ++f) {
346*9f074e33SMatthew G Knepley     const PetscInt *cone, *cellFaces;
347*9f074e33SMatthew G Knepley     PetscInt        coneSize, coff, numCellFaces, faceSize, cf;
348*9f074e33SMatthew G Knepley 
349*9f074e33SMatthew G Knepley     ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr);
350*9f074e33SMatthew G Knepley     ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr);
351*9f074e33SMatthew G Knepley     ierr = PetscSectionGetOffset(mesh->coneSection, f, &coff);CHKERRQ(ierr);
352*9f074e33SMatthew G Knepley     ierr = DMPlexGetFaces(fdm, f, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
353*9f074e33SMatthew G Knepley     if (coneSize != numCellFaces) SETERRQ3(PetscObjectComm((PetscObject)idm), PETSC_ERR_PLIB, "Invalid number of edges %D for face %D should be %D", coneSize, f, numCellFaces);
354*9f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
355*9f074e33SMatthew G Knepley       const PetscInt *econe;
356*9f074e33SMatthew G Knepley       PetscInt        esize;
357*9f074e33SMatthew G Knepley 
358*9f074e33SMatthew G Knepley       ierr = DMPlexGetConeSize(idm, cone[cf], &esize);CHKERRQ(ierr);
359*9f074e33SMatthew G Knepley       ierr = DMPlexGetCone(idm, cone[cf], &econe);CHKERRQ(ierr);
360*9f074e33SMatthew G Knepley       if (esize != 2) SETERRQ2(PetscObjectComm((PetscObject)idm), PETSC_ERR_PLIB, "Invalid number of edge endpoints %D for edge %D should be 2", esize, cone[cf]);
361*9f074e33SMatthew G Knepley       if ((cellFaces[cf*faceSize+0] == econe[0]) && (cellFaces[cf*faceSize+1] == econe[1])) {
362*9f074e33SMatthew G Knepley         /* Correctly oriented */
363*9f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = 0;
364*9f074e33SMatthew G Knepley       } else if ((cellFaces[cf*faceSize+0] == econe[1]) && (cellFaces[cf*faceSize+1] == econe[0])) {
365*9f074e33SMatthew G Knepley         /* Start at index 1, and reverse orientation */
366*9f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = -(1+1);
367*9f074e33SMatthew G Knepley       }
368*9f074e33SMatthew G Knepley     }
369*9f074e33SMatthew G Knepley   }
370*9f074e33SMatthew G Knepley   ierr = DMDestroy(&fdm);CHKERRQ(ierr);
371*9f074e33SMatthew G Knepley   /* Orient faces */
372*9f074e33SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
373*9f074e33SMatthew G Knepley     const PetscInt *cone, *cellFaces;
374*9f074e33SMatthew G Knepley     PetscInt        coneSize, coff, numCellFaces, faceSize, cf;
375*9f074e33SMatthew G Knepley 
376*9f074e33SMatthew G Knepley     ierr = DMPlexGetConeSize(idm, c, &coneSize);CHKERRQ(ierr);
377*9f074e33SMatthew G Knepley     ierr = DMPlexGetCone(idm, c, &cone);CHKERRQ(ierr);
378*9f074e33SMatthew G Knepley     ierr = PetscSectionGetOffset(mesh->coneSection, c, &coff);CHKERRQ(ierr);
379*9f074e33SMatthew G Knepley     ierr = DMPlexGetFaces(dm, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
380*9f074e33SMatthew G Knepley     if (coneSize != numCellFaces) SETERRQ3(PetscObjectComm((PetscObject)idm), PETSC_ERR_PLIB, "Invalid number of edges %D for cell %D should be %D", coneSize, c, numCellFaces);
381*9f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
382*9f074e33SMatthew G Knepley       PetscInt *origClosure = NULL, *closure;
383*9f074e33SMatthew G Knepley       PetscInt closureSize, i;
384*9f074e33SMatthew G Knepley 
385*9f074e33SMatthew G Knepley       ierr = DMPlexGetTransitiveClosure(idm, cone[cf], PETSC_TRUE, &closureSize, &origClosure);CHKERRQ(ierr);
386*9f074e33SMatthew G Knepley       if (closureSize != 7) SETERRQ2(PetscObjectComm((PetscObject)idm), PETSC_ERR_PLIB, "Invalid closure size %D for face %D should be 7", closureSize, cone[cf]);
387*9f074e33SMatthew G Knepley       for (i = 4; i < 7; ++i) {
388*9f074e33SMatthew G Knepley         if ((origClosure[i*2] < vStart) || (origClosure[i*2] >= vEnd)) SETERRQ3(PetscObjectComm((PetscObject)idm), PETSC_ERR_PLIB, "Invalid closure point %D should be a vertex in [%D, %D)", origClosure[i*2], vStart, vEnd);
389*9f074e33SMatthew G Knepley       }
390*9f074e33SMatthew G Knepley       closure = &origClosure[4*2];
391*9f074e33SMatthew G Knepley       /* Remember that this is the orientation for edges, not vertices */
392*9f074e33SMatthew G Knepley       if        ((cellFaces[cf*faceSize+0] == closure[0*2]) && (cellFaces[cf*faceSize+1] == closure[1*2]) && (cellFaces[cf*faceSize+2] == closure[2*2])) {
393*9f074e33SMatthew G Knepley         /* Correctly oriented */
394*9f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = 0;
395*9f074e33SMatthew G Knepley       } else if ((cellFaces[cf*faceSize+0] == closure[1*2]) && (cellFaces[cf*faceSize+1] == closure[2*2]) && (cellFaces[cf*faceSize+2] == closure[0*2])) {
396*9f074e33SMatthew G Knepley         /* Shifted by 1 */
397*9f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = 1;
398*9f074e33SMatthew G Knepley       } else if ((cellFaces[cf*faceSize+0] == closure[2*2]) && (cellFaces[cf*faceSize+1] == closure[0*2]) && (cellFaces[cf*faceSize+2] == closure[1*2])) {
399*9f074e33SMatthew G Knepley         /* Shifted by 2 */
400*9f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = 2;
401*9f074e33SMatthew G Knepley       } else if ((cellFaces[cf*faceSize+0] == closure[2*2]) && (cellFaces[cf*faceSize+1] == closure[1*2]) && (cellFaces[cf*faceSize+2] == closure[0*2])) {
402*9f074e33SMatthew G Knepley         /* Start at edge 1, and reverse orientation */
403*9f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = -(1+1);
404*9f074e33SMatthew G Knepley       } else if ((cellFaces[cf*faceSize+0] == closure[1*2]) && (cellFaces[cf*faceSize+1] == closure[0*2]) && (cellFaces[cf*faceSize+2] == closure[2*2])) {
405*9f074e33SMatthew G Knepley         /* Start at index 0, and reverse orientation */
406*9f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = -(0+1);
407*9f074e33SMatthew G Knepley       } else if ((cellFaces[cf*faceSize+0] == closure[0*2]) && (cellFaces[cf*faceSize+1] == closure[2*2]) && (cellFaces[cf*faceSize+2] == closure[1*2])) {
408*9f074e33SMatthew G Knepley         /* Start at index 2, and reverse orientation */
409*9f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = -(2+1);
410*9f074e33SMatthew G Knepley       } else SETERRQ3(PetscObjectComm((PetscObject)idm), PETSC_ERR_PLIB, "Face %D did not match local face %D in cell %D for any orientation", cone[cf], cf, c);
411*9f074e33SMatthew G Knepley       ierr = DMPlexRestoreTransitiveClosure(idm, cone[cf], PETSC_TRUE, &closureSize, &origClosure);CHKERRQ(ierr);
412*9f074e33SMatthew G Knepley     }
413*9f074e33SMatthew G Knepley   }
414*9f074e33SMatthew G Knepley   {
415*9f074e33SMatthew G Knepley     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr);
416*9f074e33SMatthew G Knepley     ierr = DMView(idm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
417*9f074e33SMatthew G Knepley     ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
418*9f074e33SMatthew G Knepley   }
419*9f074e33SMatthew G Knepley   *dmInt = idm;
420*9f074e33SMatthew G Knepley   PetscFunctionReturn(0);
421*9f074e33SMatthew G Knepley }
422*9f074e33SMatthew G Knepley 
423*9f074e33SMatthew G Knepley #undef __FUNCT__
424*9f074e33SMatthew G Knepley #define __FUNCT__ "DMPlexInterpolate"
425*9f074e33SMatthew G Knepley PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt)
426*9f074e33SMatthew G Knepley {
427*9f074e33SMatthew G Knepley   PetscInt       dim;
428*9f074e33SMatthew G Knepley   PetscErrorCode ierr;
429*9f074e33SMatthew G Knepley 
430*9f074e33SMatthew G Knepley   PetscFunctionBegin;
431*9f074e33SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
432*9f074e33SMatthew G Knepley   switch (dim) {
433*9f074e33SMatthew G Knepley   case 2:
434*9f074e33SMatthew G Knepley     ierr = DMPlexInterpolate_2D(dm, dmInt);CHKERRQ(ierr);break;
435*9f074e33SMatthew G Knepley   case 3:
436*9f074e33SMatthew G Knepley     ierr = DMPlexInterpolate_3D(dm, dmInt);CHKERRQ(ierr);break;
437*9f074e33SMatthew G Knepley   default:
438*9f074e33SMatthew G Knepley     SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No mesh interpolation support for dimension %D", dim);
439*9f074e33SMatthew G Knepley   }
440*9f074e33SMatthew G Knepley   PetscFunctionReturn(0);
441*9f074e33SMatthew G Knepley }
442