xref: /petsc/src/dm/impls/plex/plexinterpolate.c (revision 99836e0e8341fbabb66d730be3832b655bd5cee1)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2e8f14785SLisandro Dalcin #include <petsc/private/hashmapi.h>
3e8f14785SLisandro Dalcin #include <petsc/private/hashmapij.h>
4e8f14785SLisandro Dalcin 
5e8f14785SLisandro Dalcin /* HashIJKL */
6e8f14785SLisandro Dalcin 
7e8f14785SLisandro Dalcin #include <petsc/private/hashmap.h>
8e8f14785SLisandro Dalcin 
9e8f14785SLisandro Dalcin typedef struct _PetscHashIJKLKey { PetscInt i, j, k, l; } PetscHashIJKLKey;
10e8f14785SLisandro Dalcin 
11e8f14785SLisandro Dalcin #define PetscHashIJKLKeyHash(key) \
12e8f14785SLisandro Dalcin   PetscHashCombine(PetscHashCombine(PetscHashInt((key).i),PetscHashInt((key).j)), \
13e8f14785SLisandro Dalcin                    PetscHashCombine(PetscHashInt((key).k),PetscHashInt((key).l)))
14e8f14785SLisandro Dalcin 
15e8f14785SLisandro Dalcin #define PetscHashIJKLKeyEqual(k1,k2) \
16e8f14785SLisandro Dalcin   (((k1).i==(k2).i) ? ((k1).j==(k2).j) ? ((k1).k==(k2).k) ? ((k1).l==(k2).l) : 0 : 0 : 0)
17e8f14785SLisandro Dalcin 
18e8f14785SLisandro Dalcin PETSC_HASH_MAP(HashIJKL, PetscHashIJKLKey, PetscInt, PetscHashIJKLKeyHash, PetscHashIJKLKeyEqual, -1)
19e8f14785SLisandro Dalcin 
209f074e33SMatthew G Knepley 
219f074e33SMatthew G Knepley /*
229a5b13a2SMatthew G Knepley   DMPlexGetFaces_Internal - Gets groups of vertices that correspond to faces for the given cell
23cd921712SStefano Zampini   This assumes that the mesh is not interpolated from the depth of point p to the vertices
249f074e33SMatthew G Knepley */
25439ece16SMatthew G. Knepley PetscErrorCode DMPlexGetFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[])
269f074e33SMatthew G Knepley {
279f074e33SMatthew G Knepley   const PetscInt *cone = NULL;
28cd921712SStefano Zampini   PetscInt        coneSize;
299f074e33SMatthew G Knepley   PetscErrorCode  ierr;
309f074e33SMatthew G Knepley 
319f074e33SMatthew G Knepley   PetscFunctionBegin;
329f074e33SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
339f074e33SMatthew G Knepley   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
349f074e33SMatthew G Knepley   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
35439ece16SMatthew G. Knepley   ierr = DMPlexGetRawFaces_Internal(dm, dim, coneSize, cone, numFaces, faceSize, faces);CHKERRQ(ierr);
36439ece16SMatthew G. Knepley   PetscFunctionReturn(0);
37439ece16SMatthew G. Knepley }
38439ece16SMatthew G. Knepley 
39439ece16SMatthew G. Knepley /*
40439ece16SMatthew G. Knepley   DMPlexRestoreFaces_Internal - Restores the array
41439ece16SMatthew G. Knepley */
42439ece16SMatthew G. Knepley PetscErrorCode DMPlexRestoreFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[])
43439ece16SMatthew G. Knepley {
44439ece16SMatthew G. Knepley   PetscErrorCode  ierr;
45439ece16SMatthew G. Knepley 
46439ece16SMatthew G. Knepley   PetscFunctionBegin;
47cd921712SStefano Zampini   if (faces) { ierr = DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faces);CHKERRQ(ierr); }
48439ece16SMatthew G. Knepley   PetscFunctionReturn(0);
49439ece16SMatthew G. Knepley }
50439ece16SMatthew G. Knepley 
51439ece16SMatthew G. Knepley /*
52439ece16SMatthew G. Knepley   DMPlexGetRawFaces_Internal - Gets groups of vertices that correspond to faces for the given cone
53439ece16SMatthew G. Knepley */
54439ece16SMatthew G. Knepley PetscErrorCode DMPlexGetRawFaces_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[])
55439ece16SMatthew G. Knepley {
56439ece16SMatthew G. Knepley   PetscInt       *facesTmp;
57439ece16SMatthew G. Knepley   PetscInt        maxConeSize, maxSupportSize;
58439ece16SMatthew G. Knepley   PetscErrorCode  ierr;
59439ece16SMatthew G. Knepley 
60439ece16SMatthew G. Knepley   PetscFunctionBegin;
61439ece16SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
62cd921712SStefano Zampini   if (faces && coneSize) PetscValidIntPointer(cone,4);
63439ece16SMatthew G. Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
6469291d52SBarry Smith   if (faces) {ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), MPIU_INT, &facesTmp);CHKERRQ(ierr);}
659f074e33SMatthew G Knepley   switch (dim) {
66c49d9212SMatthew G. Knepley   case 1:
67c49d9212SMatthew G. Knepley     switch (coneSize) {
68c49d9212SMatthew G. Knepley     case 2:
69c49d9212SMatthew G. Knepley       if (faces) {
70c49d9212SMatthew G. Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
71c49d9212SMatthew G. Knepley         *faces = facesTmp;
72c49d9212SMatthew G. Knepley       }
73c49d9212SMatthew G. Knepley       if (numFaces) *numFaces = 2;
74c49d9212SMatthew G. Knepley       if (faceSize) *faceSize = 1;
75c49d9212SMatthew G. Knepley       break;
76c49d9212SMatthew G. Knepley     default:
77*99836e0eSStefano Zampini       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
78c49d9212SMatthew G. Knepley     }
79c49d9212SMatthew G. Knepley     break;
809f074e33SMatthew G Knepley   case 2:
819f074e33SMatthew G Knepley     switch (coneSize) {
829f074e33SMatthew G Knepley     case 3:
839a5b13a2SMatthew G Knepley       if (faces) {
849a5b13a2SMatthew G Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
859a5b13a2SMatthew G Knepley         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
869a5b13a2SMatthew G Knepley         facesTmp[4] = cone[2]; facesTmp[5] = cone[0];
879a5b13a2SMatthew G Knepley         *faces = facesTmp;
889a5b13a2SMatthew G Knepley       }
899a5b13a2SMatthew G Knepley       if (numFaces) *numFaces = 3;
909a5b13a2SMatthew G Knepley       if (faceSize) *faceSize = 2;
919f074e33SMatthew G Knepley       break;
929f074e33SMatthew G Knepley     case 4:
939a5b13a2SMatthew G Knepley       /* Vertices follow right hand rule */
949a5b13a2SMatthew G Knepley       if (faces) {
959a5b13a2SMatthew G Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
969a5b13a2SMatthew G Knepley         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
979a5b13a2SMatthew G Knepley         facesTmp[4] = cone[2]; facesTmp[5] = cone[3];
989a5b13a2SMatthew G Knepley         facesTmp[6] = cone[3]; facesTmp[7] = cone[0];
999a5b13a2SMatthew G Knepley         *faces = facesTmp;
1009a5b13a2SMatthew G Knepley       }
1019a5b13a2SMatthew G Knepley       if (numFaces) *numFaces = 4;
1029a5b13a2SMatthew G Knepley       if (faceSize) *faceSize = 2;
1039f074e33SMatthew G Knepley       break;
1049f074e33SMatthew G Knepley     default:
105*99836e0eSStefano Zampini       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
1069f074e33SMatthew G Knepley     }
1079f074e33SMatthew G Knepley     break;
1089f074e33SMatthew G Knepley   case 3:
1099f074e33SMatthew G Knepley     switch (coneSize) {
1109f074e33SMatthew G Knepley     case 3:
1119a5b13a2SMatthew G Knepley       if (faces) {
1129a5b13a2SMatthew G Knepley         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
1139a5b13a2SMatthew G Knepley         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
1149a5b13a2SMatthew G Knepley         facesTmp[4] = cone[2]; facesTmp[5] = cone[0];
1159a5b13a2SMatthew G Knepley         *faces = facesTmp;
1169a5b13a2SMatthew G Knepley       }
1179a5b13a2SMatthew G Knepley       if (numFaces) *numFaces = 3;
1189a5b13a2SMatthew G Knepley       if (faceSize) *faceSize = 2;
1199f074e33SMatthew G Knepley       break;
1209f074e33SMatthew G Knepley     case 4:
1212e1b13c2SMatthew G. Knepley       /* Vertices of first face follow right hand rule and normal points away from last vertex */
1229a5b13a2SMatthew G Knepley       if (faces) {
1232e1b13c2SMatthew G. Knepley         facesTmp[0] = cone[0]; facesTmp[1]  = cone[1]; facesTmp[2]  = cone[2];
1242e1b13c2SMatthew G. Knepley         facesTmp[3] = cone[0]; facesTmp[4]  = cone[3]; facesTmp[5]  = cone[1];
1252e1b13c2SMatthew G. Knepley         facesTmp[6] = cone[0]; facesTmp[7]  = cone[2]; facesTmp[8]  = cone[3];
1262e1b13c2SMatthew G. Knepley         facesTmp[9] = cone[2]; facesTmp[10] = cone[1]; facesTmp[11] = cone[3];
1279a5b13a2SMatthew G Knepley         *faces = facesTmp;
1289a5b13a2SMatthew G Knepley       }
1299a5b13a2SMatthew G Knepley       if (numFaces) *numFaces = 4;
1309a5b13a2SMatthew G Knepley       if (faceSize) *faceSize = 3;
1319a5b13a2SMatthew G Knepley       break;
1329a5b13a2SMatthew G Knepley     case 8:
133e368ef20SMatthew G. Knepley       /*  7--------6
134e368ef20SMatthew G. Knepley          /|       /|
135e368ef20SMatthew G. Knepley         / |      / |
136e368ef20SMatthew G. Knepley        4--------5  |
137e368ef20SMatthew G. Knepley        |  |     |  |
138e368ef20SMatthew G. Knepley        |  |     |  |
139e368ef20SMatthew G. Knepley        |  1--------2
140e368ef20SMatthew G. Knepley        | /      | /
141e368ef20SMatthew G. Knepley        |/       |/
142e368ef20SMatthew G. Knepley        0--------3
143e368ef20SMatthew G. Knepley        */
1449a5b13a2SMatthew G Knepley       if (faces) {
14551cfd6a4SMatthew G. Knepley         facesTmp[0]  = cone[0]; facesTmp[1]  = cone[1]; facesTmp[2]  = cone[2]; facesTmp[3]  = cone[3]; /* Bottom */
14651cfd6a4SMatthew G. Knepley         facesTmp[4]  = cone[4]; facesTmp[5]  = cone[5]; facesTmp[6]  = cone[6]; facesTmp[7]  = cone[7]; /* Top */
14751cfd6a4SMatthew G. Knepley         facesTmp[8]  = cone[0]; facesTmp[9]  = cone[3]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4]; /* Front */
14851cfd6a4SMatthew G. Knepley         facesTmp[12] = cone[2]; facesTmp[13] = cone[1]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6]; /* Back */
14951cfd6a4SMatthew G. Knepley         facesTmp[16] = cone[3]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5]; /* Right */
15051cfd6a4SMatthew G. Knepley         facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[1]; /* Left */
1519a5b13a2SMatthew G Knepley         *faces = facesTmp;
1529a5b13a2SMatthew G Knepley       }
1539a5b13a2SMatthew G Knepley       if (numFaces) *numFaces = 6;
1549a5b13a2SMatthew G Knepley       if (faceSize) *faceSize = 4;
1559f074e33SMatthew G Knepley       break;
1569f074e33SMatthew G Knepley     default:
157*99836e0eSStefano Zampini       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
1589f074e33SMatthew G Knepley     }
1599f074e33SMatthew G Knepley     break;
1609f074e33SMatthew G Knepley   default:
161*99836e0eSStefano Zampini     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim);
1629f074e33SMatthew G Knepley   }
1639f074e33SMatthew G Knepley   PetscFunctionReturn(0);
1649f074e33SMatthew G Knepley }
1659f074e33SMatthew G Knepley 
166*99836e0eSStefano Zampini /*
167*99836e0eSStefano Zampini   DMPlexGetRawFacesHybrid_Internal - Gets groups of vertices that correspond to faces for the given cone using hybrid ordering (prisms)
168*99836e0eSStefano Zampini */
169*99836e0eSStefano Zampini static PetscErrorCode DMPlexGetRawFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[])
170*99836e0eSStefano Zampini {
171*99836e0eSStefano Zampini   PetscInt       *facesTmp;
172*99836e0eSStefano Zampini   PetscInt        maxConeSize, maxSupportSize;
173*99836e0eSStefano Zampini   PetscErrorCode  ierr;
174*99836e0eSStefano Zampini 
175*99836e0eSStefano Zampini   PetscFunctionBegin;
176*99836e0eSStefano Zampini   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
177*99836e0eSStefano Zampini   if (faces && coneSize) PetscValidIntPointer(cone,4);
178*99836e0eSStefano Zampini   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
179*99836e0eSStefano Zampini   if (faces) {ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), MPIU_INT, &facesTmp);CHKERRQ(ierr);}
180*99836e0eSStefano Zampini   switch (dim) {
181*99836e0eSStefano Zampini   case 1:
182*99836e0eSStefano Zampini     switch (coneSize) {
183*99836e0eSStefano Zampini     case 2:
184*99836e0eSStefano Zampini       if (faces) {
185*99836e0eSStefano Zampini         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
186*99836e0eSStefano Zampini         *faces = facesTmp;
187*99836e0eSStefano Zampini       }
188*99836e0eSStefano Zampini       if (numFaces)     *numFaces = 2;
189*99836e0eSStefano Zampini       if (numFacesNotH) *numFacesNotH = 2;
190*99836e0eSStefano Zampini       if (faceSize)     *faceSize = 1;
191*99836e0eSStefano Zampini       break;
192*99836e0eSStefano Zampini     default:
193*99836e0eSStefano Zampini       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
194*99836e0eSStefano Zampini     }
195*99836e0eSStefano Zampini     break;
196*99836e0eSStefano Zampini   case 2:
197*99836e0eSStefano Zampini     switch (coneSize) {
198*99836e0eSStefano Zampini     case 4:
199*99836e0eSStefano Zampini       if (faces) {
200*99836e0eSStefano Zampini         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
201*99836e0eSStefano Zampini         facesTmp[2] = cone[2]; facesTmp[3] = cone[3];
202*99836e0eSStefano Zampini         facesTmp[4] = cone[0]; facesTmp[5] = cone[2];
203*99836e0eSStefano Zampini         facesTmp[6] = cone[1]; facesTmp[7] = cone[3];
204*99836e0eSStefano Zampini         *faces = facesTmp;
205*99836e0eSStefano Zampini       }
206*99836e0eSStefano Zampini       if (numFaces)     *numFaces = 4;
207*99836e0eSStefano Zampini       if (numFacesNotH) *numFacesNotH = 2;
208*99836e0eSStefano Zampini       if (faceSize)     *faceSize = 2;
209*99836e0eSStefano Zampini       break;
210*99836e0eSStefano Zampini     default:
211*99836e0eSStefano Zampini       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
212*99836e0eSStefano Zampini     }
213*99836e0eSStefano Zampini     break;
214*99836e0eSStefano Zampini   case 3:
215*99836e0eSStefano Zampini     switch (coneSize) {
216*99836e0eSStefano Zampini     case 6: /* triangular prism */
217*99836e0eSStefano Zampini       if (faces) {
218*99836e0eSStefano Zampini         facesTmp[0]  = cone[0]; facesTmp[1]  = cone[1]; facesTmp[2]  = cone[2]; facesTmp[3]  = -1;      /* Bottom */
219*99836e0eSStefano Zampini         facesTmp[4]  = cone[3]; facesTmp[5]  = cone[4]; facesTmp[6]  = cone[5]; facesTmp[7]  = -1;      /* Top */
220*99836e0eSStefano Zampini         facesTmp[8]  = cone[0]; facesTmp[9]  = cone[1]; facesTmp[10] = cone[3]; facesTmp[11] = cone[4]; /* Back left */
221*99836e0eSStefano Zampini         facesTmp[12] = cone[1]; facesTmp[13] = cone[2]; facesTmp[14] = cone[4]; facesTmp[15] = cone[5]; /* Back right */
222*99836e0eSStefano Zampini         facesTmp[16] = cone[2]; facesTmp[17] = cone[0]; facesTmp[18] = cone[5]; facesTmp[19] = cone[3]; /* Front */
223*99836e0eSStefano Zampini         *faces = facesTmp;
224*99836e0eSStefano Zampini       }
225*99836e0eSStefano Zampini       if (numFaces)     *numFaces = 5;
226*99836e0eSStefano Zampini       if (numFacesNotH) *numFacesNotH = 2;
227*99836e0eSStefano Zampini       if (faceSize)     *faceSize = -4;
228*99836e0eSStefano Zampini       break;
229*99836e0eSStefano Zampini     default:
230*99836e0eSStefano Zampini       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
231*99836e0eSStefano Zampini     }
232*99836e0eSStefano Zampini     break;
233*99836e0eSStefano Zampini   default:
234*99836e0eSStefano Zampini     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim);
235*99836e0eSStefano Zampini   }
236*99836e0eSStefano Zampini   PetscFunctionReturn(0);
237*99836e0eSStefano Zampini }
238*99836e0eSStefano Zampini 
239*99836e0eSStefano Zampini static PetscErrorCode DMPlexRestoreRawFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[])
240*99836e0eSStefano Zampini {
241*99836e0eSStefano Zampini   PetscErrorCode  ierr;
242*99836e0eSStefano Zampini 
243*99836e0eSStefano Zampini   PetscFunctionBegin;
244*99836e0eSStefano Zampini   if (faces) { ierr = DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faces);CHKERRQ(ierr); }
245*99836e0eSStefano Zampini   PetscFunctionReturn(0);
246*99836e0eSStefano Zampini }
247*99836e0eSStefano Zampini 
248*99836e0eSStefano Zampini static PetscErrorCode DMPlexGetFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[])
249*99836e0eSStefano Zampini {
250*99836e0eSStefano Zampini   const PetscInt *cone = NULL;
251*99836e0eSStefano Zampini   PetscInt        coneSize;
252*99836e0eSStefano Zampini   PetscErrorCode  ierr;
253*99836e0eSStefano Zampini 
254*99836e0eSStefano Zampini   PetscFunctionBegin;
255*99836e0eSStefano Zampini   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
256*99836e0eSStefano Zampini   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
257*99836e0eSStefano Zampini   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
258*99836e0eSStefano Zampini   ierr = DMPlexGetRawFacesHybrid_Internal(dm, dim, coneSize, cone, numFaces, numFacesNotH, faceSize, faces);CHKERRQ(ierr);
259*99836e0eSStefano Zampini   PetscFunctionReturn(0);
260*99836e0eSStefano Zampini }
261*99836e0eSStefano Zampini 
2629a5b13a2SMatthew G Knepley /* This interpolates faces for cells at some stratum */
2639a5b13a2SMatthew G Knepley static PetscErrorCode DMPlexInterpolateFaces_Internal(DM dm, PetscInt cellDepth, DM idm)
2649f074e33SMatthew G Knepley {
265d4efc6f1SMatthew G. Knepley   DMLabel        subpointMap;
2669a5b13a2SMatthew G Knepley   PetscHashIJKL  faceTable;
2679a5b13a2SMatthew G Knepley   PetscInt      *pStart, *pEnd;
2689a5b13a2SMatthew G Knepley   PetscInt       cellDim, depth, faceDepth = cellDepth, numPoints = 0, faceSizeAll = 0, face, c, d;
269*99836e0eSStefano Zampini   PetscInt       coneSizeH = 0, faceSizeAllH = 0, numCellFacesH = 0, faceH, pMax = -1, dim, outerloop;
270*99836e0eSStefano Zampini   PetscInt       cMax, fMax, eMax, vMax;
2719f074e33SMatthew G Knepley   PetscErrorCode ierr;
2729f074e33SMatthew G Knepley 
2739f074e33SMatthew G Knepley   PetscFunctionBegin;
274c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &cellDim);CHKERRQ(ierr);
275d4efc6f1SMatthew G. Knepley   /* HACK: I need a better way to determine face dimension, or an alternative to GetFaces() */
276d4efc6f1SMatthew G. Knepley   ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr);
277d4efc6f1SMatthew G. Knepley   if (subpointMap) ++cellDim;
2789a5b13a2SMatthew G Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2799a5b13a2SMatthew G Knepley   ++depth;
2809a5b13a2SMatthew G Knepley   ++cellDepth;
2819a5b13a2SMatthew G Knepley   cellDim -= depth - cellDepth;
282dcca6d9dSJed Brown   ierr = PetscMalloc2(depth+1,&pStart,depth+1,&pEnd);CHKERRQ(ierr);
2839a5b13a2SMatthew G Knepley   for (d = depth-1; d >= faceDepth; --d) {
2849a5b13a2SMatthew G Knepley     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d+1], &pEnd[d+1]);CHKERRQ(ierr);
2859f074e33SMatthew G Knepley   }
2869a5b13a2SMatthew G Knepley   ierr = DMPlexGetDepthStratum(dm, -1, NULL, &pStart[faceDepth]);CHKERRQ(ierr);
2879a5b13a2SMatthew G Knepley   pEnd[faceDepth] = pStart[faceDepth];
2889a5b13a2SMatthew G Knepley   for (d = faceDepth-1; d >= 0; --d) {
2899a5b13a2SMatthew G Knepley     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
2909f074e33SMatthew G Knepley   }
291*99836e0eSStefano Zampini   cMax = fMax = eMax = vMax = PETSC_DETERMINE;
292*99836e0eSStefano Zampini   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
293*99836e0eSStefano Zampini   if (cellDim == dim) {
294*99836e0eSStefano Zampini     ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
295*99836e0eSStefano Zampini     pMax = cMax;
296*99836e0eSStefano Zampini   } else if (cellDim == dim -1) {
297*99836e0eSStefano Zampini     ierr = DMPlexGetHybridBounds(dm, &cMax, &fMax, NULL, NULL);CHKERRQ(ierr);
298*99836e0eSStefano Zampini     pMax = fMax;
299*99836e0eSStefano Zampini   }
300*99836e0eSStefano Zampini   pMax = pMax < 0 ? pEnd[cellDepth] : pMax;
301*99836e0eSStefano Zampini   if (pMax < pEnd[cellDepth]) {
302*99836e0eSStefano Zampini     const PetscInt *cellFaces, *cone;
303*99836e0eSStefano Zampini     PetscInt        numCellFacesT, faceSize, cf;
3049f074e33SMatthew G Knepley 
305*99836e0eSStefano Zampini     ierr = DMPlexGetConeSize(dm, pMax, &coneSizeH);CHKERRQ(ierr);
306*99836e0eSStefano Zampini     ierr = DMPlexGetCone(dm, pMax, &cone);CHKERRQ(ierr);
307*99836e0eSStefano Zampini     ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSizeH, cone, &numCellFacesH, &numCellFacesT, &faceSize, &cellFaces);CHKERRQ(ierr);
308*99836e0eSStefano Zampini     if (faceSize < 0) {
309*99836e0eSStefano Zampini       PetscInt *sizes, minv, maxv;
310*99836e0eSStefano Zampini 
311*99836e0eSStefano Zampini       /* count vertices of hybrid and non-hybrid faces */
312*99836e0eSStefano Zampini       ierr = PetscCalloc1(numCellFacesH, &sizes);CHKERRQ(ierr);
313*99836e0eSStefano Zampini       for (cf = 0; cf < numCellFacesT; ++cf) { /* These are the non-hybrid faces */
314*99836e0eSStefano Zampini         const PetscInt *cellFace = &cellFaces[-cf*faceSize];
315*99836e0eSStefano Zampini         PetscInt       f;
316*99836e0eSStefano Zampini 
317*99836e0eSStefano Zampini         for (f = 0; f < -faceSize; ++f) sizes[cf] += (cellFace[f] >= 0 ? 1 : 0);
318*99836e0eSStefano Zampini       }
319*99836e0eSStefano Zampini       ierr = PetscSortInt(numCellFacesT, sizes);CHKERRQ(ierr);
320*99836e0eSStefano Zampini       minv = sizes[0];
321*99836e0eSStefano Zampini       maxv = sizes[PetscMax(numCellFacesT-1, 0)];
322*99836e0eSStefano Zampini       if (minv != maxv) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Different number of vertices for non-hybrid face %D != %D", minv, maxv);
323*99836e0eSStefano Zampini       faceSizeAll = minv;
324*99836e0eSStefano Zampini       ierr = PetscMemzero(sizes, numCellFacesH*sizeof(PetscInt));CHKERRQ(ierr);
325*99836e0eSStefano Zampini       for (cf = numCellFacesT; cf < numCellFacesH; ++cf) { /* These are the hybrid faces */
326*99836e0eSStefano Zampini         const PetscInt *cellFace = &cellFaces[-cf*faceSize];
327*99836e0eSStefano Zampini         PetscInt       f;
328*99836e0eSStefano Zampini 
329*99836e0eSStefano Zampini         for (f = 0; f < -faceSize; ++f) sizes[cf-numCellFacesT] += (cellFace[f] >= 0 ? 1 : 0);
330*99836e0eSStefano Zampini       }
331*99836e0eSStefano Zampini       ierr = PetscSortInt(numCellFacesH - numCellFacesT, sizes);CHKERRQ(ierr);
332*99836e0eSStefano Zampini       minv = sizes[0];
333*99836e0eSStefano Zampini       maxv = sizes[PetscMax(numCellFacesH - numCellFacesT-1, 0)];
334*99836e0eSStefano Zampini       ierr = PetscFree(sizes);CHKERRQ(ierr);
335*99836e0eSStefano Zampini       if (minv != maxv) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Different number of vertices for hybrid face %D != %D", minv, maxv);
336*99836e0eSStefano Zampini       faceSizeAllH = minv;
337*99836e0eSStefano Zampini     } else { /* the size of the faces in hybrid cells is the same */
338*99836e0eSStefano Zampini       faceSizeAll = faceSizeAllH = faceSize;
339*99836e0eSStefano Zampini     }
340*99836e0eSStefano Zampini     ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, cone, &numCellFacesH, &numCellFacesT, &faceSize, &cellFaces);CHKERRQ(ierr);
341*99836e0eSStefano Zampini   } else if (pEnd[cellDepth] > pStart[cellDepth]) {
342*99836e0eSStefano Zampini     ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr);
343*99836e0eSStefano Zampini   }
344*99836e0eSStefano Zampini   if (faceSizeAll > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll);
345*99836e0eSStefano Zampini 
346*99836e0eSStefano Zampini   /* With hybrid grids, we first iterate on hybrid cells and start numbering the non-hybrid grids
347*99836e0eSStefano Zampini      Then, faces for non-hybrid cells are numbered.
348*99836e0eSStefano Zampini      This is to guarantee consistent orientations (all 0) of all the points in the cone of the hybrid cells */
349*99836e0eSStefano Zampini   ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr);
350*99836e0eSStefano Zampini   for (outerloop = 0, face = pStart[faceDepth]; outerloop < 2; outerloop++) {
351*99836e0eSStefano Zampini     PetscInt start, end;
352*99836e0eSStefano Zampini 
353*99836e0eSStefano Zampini     start = outerloop == 0 ? pMax : pStart[cellDepth];
354*99836e0eSStefano Zampini     end = outerloop == 0 ? pEnd[cellDepth] : pMax;
355*99836e0eSStefano Zampini     for (c = start; c < end; ++c) {
356*99836e0eSStefano Zampini       const PetscInt *cellFaces;
357*99836e0eSStefano Zampini       PetscInt        numCellFaces, faceSize, faceSizeInc, cf;
358*99836e0eSStefano Zampini 
359*99836e0eSStefano Zampini       if (c < pMax) {
3609a5b13a2SMatthew G Knepley         ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
3619a5b13a2SMatthew G Knepley         if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll);
362*99836e0eSStefano Zampini       } else { /* Hybrid cell */
363*99836e0eSStefano Zampini         const PetscInt *cone;
364*99836e0eSStefano Zampini         PetscInt        numCellFacesN, coneSize;
365*99836e0eSStefano Zampini 
366*99836e0eSStefano Zampini         ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
367*99836e0eSStefano Zampini         if (coneSize != coneSizeH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid coneSize %D != %D", coneSize, coneSizeH);
368*99836e0eSStefano Zampini         ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
369*99836e0eSStefano Zampini         ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
370*99836e0eSStefano Zampini         if (numCellFaces != numCellFacesH) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected numCellFaces %D != %D for hybrid cell %D", numCellFaces, numCellFacesH, c);
371*99836e0eSStefano Zampini         faceSize = PetscMax(faceSize, -faceSize);
372*99836e0eSStefano Zampini         if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSize);
373*99836e0eSStefano Zampini         numCellFaces = numCellFacesN; /* process only non-hybrid faces */
374*99836e0eSStefano Zampini       }
375*99836e0eSStefano Zampini       faceSizeInc = faceSize;
3769f074e33SMatthew G Knepley       for (cf = 0; cf < numCellFaces; ++cf) {
377*99836e0eSStefano Zampini         const PetscInt   *cellFace = &cellFaces[cf*faceSizeInc];
378*99836e0eSStefano Zampini         PetscInt          faceSizeH = faceSize;
3799a5b13a2SMatthew G Knepley         PetscHashIJKLKey  key;
380e8f14785SLisandro Dalcin         PetscHashIter     iter;
381e8f14785SLisandro Dalcin         PetscBool         missing;
3829f074e33SMatthew G Knepley 
383*99836e0eSStefano Zampini         if (faceSizeInc == 2) {
3849a5b13a2SMatthew G Knepley           key.i = PetscMin(cellFace[0], cellFace[1]);
3859a5b13a2SMatthew G Knepley           key.j = PetscMax(cellFace[0], cellFace[1]);
386*99836e0eSStefano Zampini           key.k = PETSC_MAX_INT;
387*99836e0eSStefano Zampini           key.l = PETSC_MAX_INT;
3889a5b13a2SMatthew G Knepley         } else {
389*99836e0eSStefano Zampini           key.i = cellFace[0];
390*99836e0eSStefano Zampini           key.j = cellFace[1];
391*99836e0eSStefano Zampini           key.k = cellFace[2];
392*99836e0eSStefano Zampini           key.l = faceSize > 3 ? (cellFace[3] < 0 ? faceSizeH = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT;
393302440fdSBarry Smith           ierr  = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr);
3949f074e33SMatthew G Knepley         }
395*99836e0eSStefano Zampini         /* this check is redundant for non-hybrid meshes */
396*99836e0eSStefano Zampini         if (faceSizeH != faceSizeAll) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected number of vertices for face %D of point %D -> %D != %D", cf, c, faceSizeH, faceSizeAll);
397e8f14785SLisandro Dalcin         ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr);
398e8f14785SLisandro Dalcin         if (missing) {ierr = PetscHashIJKLIterSet(faceTable, iter, face++);CHKERRQ(ierr);}
3999a5b13a2SMatthew G Knepley       }
400*99836e0eSStefano Zampini       if (c < pMax) {
401439ece16SMatthew G. Knepley         ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
402*99836e0eSStefano Zampini       } else {
403*99836e0eSStefano Zampini         ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, NULL, NULL, NULL, NULL, &cellFaces);CHKERRQ(ierr);
404*99836e0eSStefano Zampini       }
405*99836e0eSStefano Zampini     }
406*99836e0eSStefano Zampini   }
407*99836e0eSStefano Zampini   pEnd[faceDepth] = face;
408*99836e0eSStefano Zampini 
409*99836e0eSStefano Zampini   /* Second pass for hybrid meshes: number hybrid faces */
410*99836e0eSStefano Zampini   for (c = pMax; c < pEnd[cellDepth]; ++c) {
411*99836e0eSStefano Zampini     const PetscInt *cellFaces, *cone;
412*99836e0eSStefano Zampini     PetscInt        numCellFaces, numCellFacesN, faceSize, cf, coneSize;
413*99836e0eSStefano Zampini 
414*99836e0eSStefano Zampini     ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
415*99836e0eSStefano Zampini     ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
416*99836e0eSStefano Zampini     ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
417*99836e0eSStefano Zampini     if (numCellFaces != numCellFacesH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid numCellFaces %D != %D", numCellFaces, numCellFacesH);
418*99836e0eSStefano Zampini     faceSize = PetscMax(faceSize, -faceSize);
419*99836e0eSStefano Zampini     for (cf = numCellFacesN; cf < numCellFaces; ++cf) { /* These are the hybrid faces */
420*99836e0eSStefano Zampini       const PetscInt   *cellFace = &cellFaces[cf*faceSize];
421*99836e0eSStefano Zampini       PetscHashIJKLKey  key;
422*99836e0eSStefano Zampini       PetscHashIter     iter;
423*99836e0eSStefano Zampini       PetscBool         missing;
424*99836e0eSStefano Zampini       PetscInt          faceSizeH = faceSize;
425*99836e0eSStefano Zampini 
426*99836e0eSStefano Zampini       if (faceSize == 2) {
427*99836e0eSStefano Zampini         key.i = PetscMin(cellFace[0], cellFace[1]);
428*99836e0eSStefano Zampini         key.j = PetscMax(cellFace[0], cellFace[1]);
429*99836e0eSStefano Zampini         key.k = PETSC_MAX_INT;
430*99836e0eSStefano Zampini         key.l = PETSC_MAX_INT;
431*99836e0eSStefano Zampini       } else {
432*99836e0eSStefano Zampini         key.i = cellFace[0];
433*99836e0eSStefano Zampini         key.j = cellFace[1];
434*99836e0eSStefano Zampini         key.k = cellFace[2];
435*99836e0eSStefano Zampini         key.l = faceSize > 3 ? (cellFace[3] < 0 ? faceSizeH = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT;
436*99836e0eSStefano Zampini         ierr  = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr);
437*99836e0eSStefano Zampini       }
438*99836e0eSStefano Zampini       if (faceSizeH != faceSizeAllH) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected number of vertices for hybrid face %D of point %D -> %D != %D", cf, c, faceSizeH, faceSizeAllH);
439*99836e0eSStefano Zampini       ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr);
440*99836e0eSStefano Zampini       if (missing) {ierr = PetscHashIJKLIterSet(faceTable, iter, face++);CHKERRQ(ierr);}
441*99836e0eSStefano Zampini     }
442*99836e0eSStefano Zampini     ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
443*99836e0eSStefano Zampini   }
444*99836e0eSStefano Zampini   faceH = face - pEnd[faceDepth];
445*99836e0eSStefano Zampini   if (faceH) {
446*99836e0eSStefano Zampini     if (fMax == PETSC_DETERMINE) fMax = pEnd[faceDepth];
447*99836e0eSStefano Zampini     else if (eMax == PETSC_DETERMINE) eMax = pEnd[faceDepth];
448*99836e0eSStefano Zampini     else SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of unassigned hybrid facets %D for cellDim %D and dimension %D", faceH, cellDim, dim);
4499a5b13a2SMatthew G Knepley   }
4509a5b13a2SMatthew G Knepley   pEnd[faceDepth] = face;
4519a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr);
4529a5b13a2SMatthew G Knepley   /* Count new points */
4539a5b13a2SMatthew G Knepley   for (d = 0; d <= depth; ++d) {
4549a5b13a2SMatthew G Knepley     numPoints += pEnd[d]-pStart[d];
4559a5b13a2SMatthew G Knepley   }
4569a5b13a2SMatthew G Knepley   ierr = DMPlexSetChart(idm, 0, numPoints);CHKERRQ(ierr);
4579a5b13a2SMatthew G Knepley   /* Set cone sizes */
4589a5b13a2SMatthew G Knepley   for (d = 0; d <= depth; ++d) {
4599a5b13a2SMatthew G Knepley     PetscInt coneSize, p;
4609f074e33SMatthew G Knepley 
4619a5b13a2SMatthew G Knepley     if (d == faceDepth) {
4629a5b13a2SMatthew G Knepley       /* I see no way to do this if we admit faces of different shapes */
463*99836e0eSStefano Zampini       for (p = pStart[d]; p < pEnd[d]-faceH; ++p) {
4649a5b13a2SMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr);
4659a5b13a2SMatthew G Knepley       }
466*99836e0eSStefano Zampini       for (p = pEnd[d]-faceH; p < pEnd[d]; ++p) {
467*99836e0eSStefano Zampini         ierr = DMPlexSetConeSize(idm, p, faceSizeAllH);CHKERRQ(ierr);
468*99836e0eSStefano Zampini       }
469a014044eSMatthew G Knepley     } else if (d == cellDepth) {
470a014044eSMatthew G Knepley       for (p = pStart[d]; p < pEnd[d]; ++p) {
471a014044eSMatthew G Knepley         /* Number of cell faces may be different from number of cell vertices*/
472*99836e0eSStefano Zampini         if (p < pMax) {
473a014044eSMatthew G Knepley           ierr = DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);CHKERRQ(ierr);
474*99836e0eSStefano Zampini         } else {
475*99836e0eSStefano Zampini           ierr = DMPlexGetFacesHybrid_Internal(dm, cellDim, p, &coneSize, NULL, NULL, NULL);CHKERRQ(ierr);
476*99836e0eSStefano Zampini         }
477a014044eSMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr);
478a014044eSMatthew G Knepley       }
4799a5b13a2SMatthew G Knepley     } else {
4809a5b13a2SMatthew G Knepley       for (p = pStart[d]; p < pEnd[d]; ++p) {
4819a5b13a2SMatthew G Knepley         ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
4829a5b13a2SMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr);
4839f074e33SMatthew G Knepley       }
4849f074e33SMatthew G Knepley     }
4859f074e33SMatthew G Knepley   }
4869f074e33SMatthew G Knepley   ierr = DMSetUp(idm);CHKERRQ(ierr);
4879f074e33SMatthew G Knepley   /* Get face cones from subsets of cell vertices */
488*99836e0eSStefano Zampini   if (faceSizeAll > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll);
4899a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr);
4909a5b13a2SMatthew G Knepley   for (d = depth; d > cellDepth; --d) {
4919a5b13a2SMatthew G Knepley     const PetscInt *cone;
4929a5b13a2SMatthew G Knepley     PetscInt        p;
4939a5b13a2SMatthew G Knepley 
4949a5b13a2SMatthew G Knepley     for (p = pStart[d]; p < pEnd[d]; ++p) {
4959a5b13a2SMatthew G Knepley       ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
4969a5b13a2SMatthew G Knepley       ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr);
4979a5b13a2SMatthew G Knepley       ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr);
4989a5b13a2SMatthew G Knepley       ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr);
4999f074e33SMatthew G Knepley     }
5009a5b13a2SMatthew G Knepley   }
501*99836e0eSStefano Zampini   for (outerloop = 0, face = pStart[faceDepth]; outerloop < 2; outerloop++) {
502*99836e0eSStefano Zampini     PetscInt start, end;
5039f074e33SMatthew G Knepley 
504*99836e0eSStefano Zampini     start = outerloop == 0 ? pMax : pStart[cellDepth];
505*99836e0eSStefano Zampini     end = outerloop == 0 ? pEnd[cellDepth] : pMax;
506*99836e0eSStefano Zampini     for (c = start; c < end; ++c) {
507*99836e0eSStefano Zampini       const PetscInt *cellFaces;
508*99836e0eSStefano Zampini       PetscInt        numCellFaces, faceSize, faceSizeInc, cf;
509*99836e0eSStefano Zampini 
510*99836e0eSStefano Zampini       if (c < pMax) {
5119a5b13a2SMatthew G Knepley         ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
5129a5b13a2SMatthew G Knepley         if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll);
513*99836e0eSStefano Zampini       } else {
514*99836e0eSStefano Zampini         const PetscInt *cone;
515*99836e0eSStefano Zampini         PetscInt        numCellFacesN, coneSize;
516*99836e0eSStefano Zampini 
517*99836e0eSStefano Zampini         ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
518*99836e0eSStefano Zampini         if (coneSize != coneSizeH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid coneSize %D != %D", coneSize, coneSizeH);
519*99836e0eSStefano Zampini         ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
520*99836e0eSStefano Zampini         ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
521*99836e0eSStefano Zampini         if (numCellFaces != numCellFacesH) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected numCellFaces %D != %D for hybrid cell %D", numCellFaces, numCellFacesH, c);
522*99836e0eSStefano Zampini         faceSize = PetscMax(faceSize, -faceSize);
523*99836e0eSStefano Zampini         if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSize);
524*99836e0eSStefano Zampini         numCellFaces = numCellFacesN; /* process only non-hybrid faces */
525*99836e0eSStefano Zampini       }
526*99836e0eSStefano Zampini       faceSizeInc = faceSize;
5279f074e33SMatthew G Knepley       for (cf = 0; cf < numCellFaces; ++cf) {
528*99836e0eSStefano Zampini         const PetscInt  *cellFace = &cellFaces[cf*faceSizeInc];
5299a5b13a2SMatthew G Knepley         PetscHashIJKLKey key;
530e8f14785SLisandro Dalcin         PetscHashIter    iter;
531e8f14785SLisandro Dalcin         PetscBool        missing;
5329f074e33SMatthew G Knepley 
533*99836e0eSStefano Zampini         if (faceSizeInc == 2) {
5349a5b13a2SMatthew G Knepley           key.i = PetscMin(cellFace[0], cellFace[1]);
5359a5b13a2SMatthew G Knepley           key.j = PetscMax(cellFace[0], cellFace[1]);
536*99836e0eSStefano Zampini           key.k = PETSC_MAX_INT;
537*99836e0eSStefano Zampini           key.l = PETSC_MAX_INT;
5389a5b13a2SMatthew G Knepley         } else {
539e8f14785SLisandro Dalcin           key.i = cellFace[0];
540e8f14785SLisandro Dalcin           key.j = cellFace[1];
541e8f14785SLisandro Dalcin           key.k = cellFace[2];
542*99836e0eSStefano Zampini           key.l = faceSizeInc > 3 ? (cellFace[3] < 0 ? faceSize = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT;
543*99836e0eSStefano Zampini           ierr  = PetscSortInt(faceSizeInc, (PetscInt *) &key);CHKERRQ(ierr);
5449f074e33SMatthew G Knepley         }
545e8f14785SLisandro Dalcin         ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr);
546735a0255SMatthew G. Knepley         if (missing) {
5479a5b13a2SMatthew G Knepley           ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr);
548e8f14785SLisandro Dalcin           ierr = PetscHashIJKLIterSet(faceTable, iter, face);CHKERRQ(ierr);
549735a0255SMatthew G. Knepley           ierr = DMPlexInsertCone(idm, c, cf, face++);CHKERRQ(ierr);
5509a5b13a2SMatthew G Knepley         } else {
5519a5b13a2SMatthew G Knepley           const PetscInt *cone;
552735a0255SMatthew G. Knepley           PetscInt        coneSize, ornt, i, j, f;
5539f074e33SMatthew G Knepley 
554e8f14785SLisandro Dalcin           ierr = PetscHashIJKLIterGet(faceTable, iter, &f);CHKERRQ(ierr);
5559a5b13a2SMatthew G Knepley           ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr);
5562e1b13c2SMatthew G. Knepley           /* Orient face: Do not allow reverse orientation at the first vertex */
5579f074e33SMatthew G Knepley           ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr);
5589f074e33SMatthew G Knepley           ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr);
5599a5b13a2SMatthew 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);
5609a5b13a2SMatthew G Knepley           /* - First find the initial vertex */
5619a5b13a2SMatthew G Knepley           for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break;
5629a5b13a2SMatthew G Knepley           /* - Try forward comparison */
5639a5b13a2SMatthew G Knepley           for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break;
5649a5b13a2SMatthew G Knepley           if (j == faceSize) {
5659a5b13a2SMatthew G Knepley             if ((faceSize == 2) && (i == 1)) ornt = -2;
5669a5b13a2SMatthew G Knepley             else                             ornt = i;
5679a5b13a2SMatthew G Knepley           } else {
5689a5b13a2SMatthew G Knepley             /* - Try backward comparison */
5699a5b13a2SMatthew G Knepley             for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break;
5702e1b13c2SMatthew G. Knepley             if (j == faceSize) {
5712e1b13c2SMatthew G. Knepley               if (i == 0) ornt = -faceSize;
572dc1a705cSMatthew G. Knepley               else        ornt = -i;
5732e1b13c2SMatthew G. Knepley             } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation");
5749a5b13a2SMatthew G Knepley           }
5759a5b13a2SMatthew G Knepley           ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr);
5769f074e33SMatthew G Knepley         }
5779f074e33SMatthew G Knepley       }
578*99836e0eSStefano Zampini       if (c < pMax) {
579439ece16SMatthew G. Knepley         ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
580*99836e0eSStefano Zampini       } else {
581*99836e0eSStefano Zampini         ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, NULL, NULL, NULL, NULL, &cellFaces);CHKERRQ(ierr);
5829f074e33SMatthew G Knepley       }
583*99836e0eSStefano Zampini     }
584*99836e0eSStefano Zampini   }
585*99836e0eSStefano Zampini   /* Second pass for hybrid meshes: orient hybrid faces */
586*99836e0eSStefano Zampini   for (c = pMax; c < pEnd[cellDepth]; ++c) {
587*99836e0eSStefano Zampini     const PetscInt *cellFaces, *cone;
588*99836e0eSStefano Zampini     PetscInt        numCellFaces, numCellFacesN, faceSize, cf, coneSize;
589*99836e0eSStefano Zampini 
590*99836e0eSStefano Zampini     ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
591*99836e0eSStefano Zampini     ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
592*99836e0eSStefano Zampini     ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
593*99836e0eSStefano Zampini     if (numCellFaces != numCellFacesH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid numCellFaces %D != %D", numCellFaces, numCellFacesH);
594*99836e0eSStefano Zampini     faceSize = PetscMax(faceSize, -faceSize);
595*99836e0eSStefano Zampini     for (cf = numCellFacesN; cf < numCellFaces; ++cf) { /* These are the hybrid faces */
596*99836e0eSStefano Zampini       const PetscInt   *cellFace = &cellFaces[cf*faceSize];
597*99836e0eSStefano Zampini       PetscHashIJKLKey key;
598*99836e0eSStefano Zampini       PetscHashIter    iter;
599*99836e0eSStefano Zampini       PetscBool        missing;
600*99836e0eSStefano Zampini       PetscInt         faceSizeH = faceSize;
601*99836e0eSStefano Zampini 
602*99836e0eSStefano Zampini       if (faceSize == 2) {
603*99836e0eSStefano Zampini         key.i = PetscMin(cellFace[0], cellFace[1]);
604*99836e0eSStefano Zampini         key.j = PetscMax(cellFace[0], cellFace[1]);
605*99836e0eSStefano Zampini         key.k = PETSC_MAX_INT;
606*99836e0eSStefano Zampini         key.l = PETSC_MAX_INT;
607*99836e0eSStefano Zampini       } else {
608*99836e0eSStefano Zampini         key.i = cellFace[0];
609*99836e0eSStefano Zampini         key.j = cellFace[1];
610*99836e0eSStefano Zampini         key.k = cellFace[2];
611*99836e0eSStefano Zampini         key.l = faceSize > 3 ? (cellFace[3] < 0 ? faceSizeH = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT;
612*99836e0eSStefano Zampini         ierr  = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr);
613*99836e0eSStefano Zampini       }
614*99836e0eSStefano Zampini       if (faceSizeH != faceSizeAllH) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected number of vertices for hybrid face %D of point %D -> %D != %D", cf, c, faceSizeH, faceSizeAllH);
615*99836e0eSStefano Zampini       ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr);
616*99836e0eSStefano Zampini       if (missing) {
617*99836e0eSStefano Zampini         ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr);
618*99836e0eSStefano Zampini         ierr = PetscHashIJKLIterSet(faceTable, iter, face);CHKERRQ(ierr);
619*99836e0eSStefano Zampini         ierr = DMPlexInsertCone(idm, c, cf, face++);CHKERRQ(ierr);
620*99836e0eSStefano Zampini       } else {
621*99836e0eSStefano Zampini         const PetscInt *cone;
622*99836e0eSStefano Zampini         PetscInt        coneSize, ornt, i, j, f;
623*99836e0eSStefano Zampini 
624*99836e0eSStefano Zampini         ierr = PetscHashIJKLIterGet(faceTable, iter, &f);CHKERRQ(ierr);
625*99836e0eSStefano Zampini         ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr);
626*99836e0eSStefano Zampini         /* Orient face: Do not allow reverse orientation at the first vertex */
627*99836e0eSStefano Zampini         ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr);
628*99836e0eSStefano Zampini         ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr);
629*99836e0eSStefano Zampini         if (coneSize != faceSizeH) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of face vertices %D for face %D should be %D", coneSize, f, faceSizeH);
630*99836e0eSStefano Zampini         /* - First find the initial vertex */
631*99836e0eSStefano Zampini         for (i = 0; i < faceSizeH; ++i) if (cellFace[0] == cone[i]) break;
632*99836e0eSStefano Zampini         /* - Try forward comparison */
633*99836e0eSStefano Zampini         for (j = 0; j < faceSizeH; ++j) if (cellFace[j] != cone[(i+j)%faceSizeH]) break;
634*99836e0eSStefano Zampini         if (j == faceSizeH) {
635*99836e0eSStefano Zampini           if ((faceSizeH == 2) && (i == 1)) ornt = -2;
636*99836e0eSStefano Zampini           else                             ornt = i;
637*99836e0eSStefano Zampini         } else {
638*99836e0eSStefano Zampini           /* - Try backward comparison */
639*99836e0eSStefano Zampini           for (j = 0; j < faceSizeH; ++j) if (cellFace[j] != cone[(i+faceSizeH-j)%faceSizeH]) break;
640*99836e0eSStefano Zampini           if (j == faceSizeH) {
641*99836e0eSStefano Zampini             if (i == 0) ornt = -faceSizeH;
642*99836e0eSStefano Zampini             else        ornt = -i;
643*99836e0eSStefano Zampini           } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation");
644*99836e0eSStefano Zampini         }
645*99836e0eSStefano Zampini         ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr);
646*99836e0eSStefano Zampini       }
647*99836e0eSStefano Zampini     }
648*99836e0eSStefano Zampini     ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
649*99836e0eSStefano Zampini   }
650*99836e0eSStefano Zampini   if (face != pEnd[faceDepth]) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of faces %D should be %D", face-pStart[faceDepth], pEnd[faceDepth]-pStart[faceDepth]);
651c907b753SJed Brown   ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr);
6529a5b13a2SMatthew G Knepley   ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr);
6536551a8c7SMatthew G. Knepley   ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr);
654*99836e0eSStefano Zampini   ierr = DMPlexSetHybridBounds(idm, cMax, fMax, eMax, vMax);CHKERRQ(ierr);
6559a5b13a2SMatthew G Knepley   ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr);
6569a5b13a2SMatthew G Knepley   ierr = DMPlexStratify(idm);CHKERRQ(ierr);
6579f074e33SMatthew G Knepley   PetscFunctionReturn(0);
6589f074e33SMatthew G Knepley }
6599f074e33SMatthew G Knepley 
6607bffde78SMichael Lange /* This interpolates the PointSF in parallel following local interpolation */
6617bffde78SMichael Lange static PetscErrorCode DMPlexInterpolatePointSF(DM dm, PetscSF pointSF, PetscInt depth)
6627bffde78SMichael Lange {
6639852e123SBarry Smith   PetscMPIInt        size, rank;
6647bffde78SMichael Lange   PetscInt           p, c, d, dof, offset;
665ffd6914dSSatish Balay   PetscInt           numLeaves, numRoots, candidatesSize, candidatesRemoteSize;
6667bffde78SMichael Lange   const PetscInt    *localPoints;
6677bffde78SMichael Lange   const PetscSFNode *remotePoints;
6687bffde78SMichael Lange   PetscSFNode       *candidates, *candidatesRemote, *claims;
6697bffde78SMichael Lange   PetscSection       candidateSection, candidateSectionRemote, claimSection;
670e8f14785SLisandro Dalcin   PetscHMapI         leafhash;
671e8f14785SLisandro Dalcin   PetscHMapIJ        roothash;
6727bffde78SMichael Lange   PetscHashIJKey     key;
6737bffde78SMichael Lange   PetscErrorCode     ierr;
6747bffde78SMichael Lange 
6757bffde78SMichael Lange   PetscFunctionBegin;
6769852e123SBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr);
6777bffde78SMichael Lange   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
6787bffde78SMichael Lange   ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr);
6799852e123SBarry Smith   if (size < 2 || numRoots < 0) PetscFunctionReturn(0);
68025afeb17SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr);
6817bffde78SMichael Lange   /* Build hashes of points in the SF for efficient lookup */
682e8f14785SLisandro Dalcin   ierr = PetscHMapICreate(&leafhash);CHKERRQ(ierr);
683e8f14785SLisandro Dalcin   ierr = PetscHMapIJCreate(&roothash);CHKERRQ(ierr);
6847bffde78SMichael Lange   for (p = 0; p < numLeaves; ++p) {
685e8f14785SLisandro Dalcin     ierr = PetscHMapISet(leafhash, localPoints[p], p);CHKERRQ(ierr);
686e8f14785SLisandro Dalcin     key.i = remotePoints[p].index;
687e8f14785SLisandro Dalcin     key.j = remotePoints[p].rank;
688e8f14785SLisandro Dalcin     ierr = PetscHMapIJSet(roothash, key, p);CHKERRQ(ierr);
6897bffde78SMichael Lange   }
6907bffde78SMichael Lange   /* Build a section / SFNode array of candidate points in the single-level adjacency of leaves,
6917bffde78SMichael Lange      where each candidate is defined by the root entry for the other vertex that defines the edge. */
6927bffde78SMichael Lange   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSection);CHKERRQ(ierr);
6937bffde78SMichael Lange   ierr = PetscSectionSetChart(candidateSection, 0, numRoots);CHKERRQ(ierr);
6947bffde78SMichael Lange   {
695ffd6914dSSatish Balay     PetscInt leaf, root, idx, a, *adj = NULL;
6967bffde78SMichael Lange     for (p = 0; p < numLeaves; ++p) {
6977bffde78SMichael Lange       PetscInt adjSize = PETSC_DETERMINE;
6987bffde78SMichael Lange       ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr);
6997bffde78SMichael Lange       for (a = 0; a < adjSize; ++a) {
700e8f14785SLisandro Dalcin         ierr = PetscHMapIGet(leafhash, adj[a], &leaf);CHKERRQ(ierr);
7017bffde78SMichael Lange         if (leaf >= 0) {ierr = PetscSectionAddDof(candidateSection, localPoints[p], 1);CHKERRQ(ierr);}
7027bffde78SMichael Lange       }
7037bffde78SMichael Lange     }
7047bffde78SMichael Lange     ierr = PetscSectionSetUp(candidateSection);CHKERRQ(ierr);
7057bffde78SMichael Lange     ierr = PetscSectionGetStorageSize(candidateSection, &candidatesSize);CHKERRQ(ierr);
7067bffde78SMichael Lange     ierr = PetscMalloc1(candidatesSize, &candidates);CHKERRQ(ierr);
7077bffde78SMichael Lange     for (p = 0; p < numLeaves; ++p) {
7087bffde78SMichael Lange       PetscInt adjSize = PETSC_DETERMINE;
7097bffde78SMichael Lange       ierr = PetscSectionGetOffset(candidateSection, localPoints[p], &offset);CHKERRQ(ierr);
7107bffde78SMichael Lange       ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr);
7117bffde78SMichael Lange       for (idx = 0, a = 0; a < adjSize; ++a) {
712e8f14785SLisandro Dalcin         ierr = PetscHMapIGet(leafhash, adj[a], &root);CHKERRQ(ierr);
7137bffde78SMichael Lange         if (root >= 0) candidates[offset+idx++] = remotePoints[root];
7147bffde78SMichael Lange       }
7157bffde78SMichael Lange     }
7167bffde78SMichael Lange     ierr = PetscFree(adj);CHKERRQ(ierr);
7177bffde78SMichael Lange   }
7187bffde78SMichael Lange   /* Gather candidate section / array pair into the root partition via inverse(multi(pointSF)). */
7197bffde78SMichael Lange   {
7207bffde78SMichael Lange     PetscSF   sfMulti, sfInverse, sfCandidates;
7217bffde78SMichael Lange     PetscInt *remoteOffsets;
7227bffde78SMichael Lange     ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr);
7237bffde78SMichael Lange     ierr = PetscSFCreateInverseSF(sfMulti, &sfInverse);CHKERRQ(ierr);
7247bffde78SMichael Lange     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSectionRemote);CHKERRQ(ierr);
7257bffde78SMichael Lange     ierr = PetscSFDistributeSection(sfInverse, candidateSection, &remoteOffsets, candidateSectionRemote);CHKERRQ(ierr);
7267bffde78SMichael Lange     ierr = PetscSFCreateSectionSF(sfInverse, candidateSection, remoteOffsets, candidateSectionRemote, &sfCandidates);CHKERRQ(ierr);
7277bffde78SMichael Lange     ierr = PetscSectionGetStorageSize(candidateSectionRemote, &candidatesRemoteSize);CHKERRQ(ierr);
7287bffde78SMichael Lange     ierr = PetscMalloc1(candidatesRemoteSize, &candidatesRemote);CHKERRQ(ierr);
7297bffde78SMichael Lange     ierr = PetscSFBcastBegin(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr);
7307bffde78SMichael Lange     ierr = PetscSFBcastEnd(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr);
7317bffde78SMichael Lange     ierr = PetscSFDestroy(&sfInverse);CHKERRQ(ierr);
7327bffde78SMichael Lange     ierr = PetscSFDestroy(&sfCandidates);CHKERRQ(ierr);
7337bffde78SMichael Lange     ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
7347bffde78SMichael Lange   }
7357bffde78SMichael Lange   /* Walk local roots and check for each remote candidate whether we know all required points,
7367bffde78SMichael Lange      either from owning it or having a root entry in the point SF. If we do we place a claim
7377bffde78SMichael Lange      by replacing the vertex number with our edge ID. */
7387bffde78SMichael Lange   {
7397bffde78SMichael Lange     PetscInt        idx, root, joinSize, vertices[2];
7407bffde78SMichael Lange     const PetscInt *rootdegree, *join = NULL;
7417bffde78SMichael Lange     ierr = PetscSFComputeDegreeBegin(pointSF, &rootdegree);CHKERRQ(ierr);
7427bffde78SMichael Lange     ierr = PetscSFComputeDegreeEnd(pointSF, &rootdegree);CHKERRQ(ierr);
7437bffde78SMichael Lange     /* Loop remote edge connections and put in a claim if both vertices are known */
7447bffde78SMichael Lange     for (idx = 0, p = 0; p < numRoots; ++p) {
7457bffde78SMichael Lange       for (d = 0; d < rootdegree[p]; ++d) {
7467bffde78SMichael Lange         ierr = PetscSectionGetDof(candidateSectionRemote, idx, &dof);CHKERRQ(ierr);
7477bffde78SMichael Lange         ierr = PetscSectionGetOffset(candidateSectionRemote, idx, &offset);CHKERRQ(ierr);
7487bffde78SMichael Lange         for (c = 0; c < dof; ++c) {
7497bffde78SMichael Lange           /* We own both vertices, so we claim the edge by replacing vertex with edge */
7507bffde78SMichael Lange           if (candidatesRemote[offset+c].rank == rank) {
7517bffde78SMichael Lange             vertices[0] = p; vertices[1] = candidatesRemote[offset+c].index;
7527bffde78SMichael Lange             ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
7537bffde78SMichael Lange             if (joinSize == 1) candidatesRemote[offset+c].index = join[0];
7547bffde78SMichael Lange             ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
7557bffde78SMichael Lange             continue;
7567bffde78SMichael Lange           }
7577bffde78SMichael Lange           /* If we own one vertex and share a root with the other, we claim it */
758e8f14785SLisandro Dalcin           key.i = candidatesRemote[offset+c].index;
759e8f14785SLisandro Dalcin           key.j = candidatesRemote[offset+c].rank;
760e8f14785SLisandro Dalcin           ierr = PetscHMapIJGet(roothash, key, &root);CHKERRQ(ierr);
7617bffde78SMichael Lange           if (root >= 0) {
7627bffde78SMichael Lange             vertices[0] = p; vertices[1] = localPoints[root];
7637bffde78SMichael Lange             ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
7647bffde78SMichael Lange             if (joinSize == 1) {
7657bffde78SMichael Lange               candidatesRemote[offset+c].index = join[0];
7667bffde78SMichael Lange               candidatesRemote[offset+c].rank = rank;
7677bffde78SMichael Lange             }
7687bffde78SMichael Lange             ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
7697bffde78SMichael Lange           }
7707bffde78SMichael Lange         }
7717bffde78SMichael Lange         idx++;
7727bffde78SMichael Lange       }
7737bffde78SMichael Lange     }
7747bffde78SMichael Lange   }
7757bffde78SMichael Lange   /* Push claims back to receiver via the MultiSF and derive new pointSF mapping on receiver */
7767bffde78SMichael Lange   {
7777bffde78SMichael Lange     PetscSF         sfMulti, sfClaims, sfPointNew;
778e8f14785SLisandro Dalcin     PetscHMapI      claimshash;
7797bffde78SMichael Lange     PetscInt        size, pStart, pEnd, root, joinSize, numLocalNew;
7807bffde78SMichael Lange     PetscInt       *remoteOffsets, *localPointsNew, vertices[2];
781ffd6914dSSatish Balay     const PetscInt *join = NULL;
7827bffde78SMichael Lange     PetscSFNode    *remotePointsNew;
7837bffde78SMichael Lange     ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr);
7847bffde78SMichael Lange     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &claimSection);CHKERRQ(ierr);
7857bffde78SMichael Lange     ierr = PetscSFDistributeSection(sfMulti, candidateSectionRemote, &remoteOffsets, claimSection);CHKERRQ(ierr);
7867bffde78SMichael Lange     ierr = PetscSFCreateSectionSF(sfMulti, candidateSectionRemote, remoteOffsets, claimSection, &sfClaims);CHKERRQ(ierr);
7877bffde78SMichael Lange     ierr = PetscSectionGetStorageSize(claimSection, &size);CHKERRQ(ierr);
7887bffde78SMichael Lange     ierr = PetscMalloc1(size, &claims);CHKERRQ(ierr);
7897bffde78SMichael Lange     ierr = PetscSFBcastBegin(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr);
7907bffde78SMichael Lange     ierr = PetscSFBcastEnd(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr);
7917bffde78SMichael Lange     ierr = PetscSFDestroy(&sfClaims);CHKERRQ(ierr);
7927bffde78SMichael Lange     ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
7937bffde78SMichael Lange     /* Walk the original section of local supports and add an SF entry for each updated item */
794e8f14785SLisandro Dalcin     ierr = PetscHMapICreate(&claimshash);CHKERRQ(ierr);
7957bffde78SMichael Lange     for (p = 0; p < numRoots; ++p) {
7967bffde78SMichael Lange       ierr = PetscSectionGetDof(candidateSection, p, &dof);CHKERRQ(ierr);
7977bffde78SMichael Lange       ierr = PetscSectionGetOffset(candidateSection, p, &offset);CHKERRQ(ierr);
7987bffde78SMichael Lange       for (d = 0; d < dof; ++d) {
7997bffde78SMichael Lange         if (candidates[offset+d].index != claims[offset+d].index) {
800e8f14785SLisandro Dalcin           key.i = candidates[offset+d].index;
801e8f14785SLisandro Dalcin           key.j = candidates[offset+d].rank;
802e8f14785SLisandro Dalcin           ierr = PetscHMapIJGet(roothash, key, &root);CHKERRQ(ierr);
8037bffde78SMichael Lange           if (root >= 0) {
8047bffde78SMichael Lange             vertices[0] = p; vertices[1] = localPoints[root];
8057bffde78SMichael Lange             ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
806e8f14785SLisandro Dalcin             if (joinSize == 1) {ierr = PetscHMapISet(claimshash, join[0], offset+d);CHKERRQ(ierr);}
8077bffde78SMichael Lange             ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
8087bffde78SMichael Lange           }
8097bffde78SMichael Lange         }
8107bffde78SMichael Lange       }
8117bffde78SMichael Lange     }
8127bffde78SMichael Lange     /* Create new pointSF from hashed claims */
813e8f14785SLisandro Dalcin     ierr = PetscHMapIGetSize(claimshash, &numLocalNew);CHKERRQ(ierr);
8147bffde78SMichael Lange     ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
8157bffde78SMichael Lange     ierr = PetscMalloc1(numLeaves + numLocalNew, &localPointsNew);CHKERRQ(ierr);
8167bffde78SMichael Lange     ierr = PetscMalloc1(numLeaves + numLocalNew, &remotePointsNew);CHKERRQ(ierr);
8177bffde78SMichael Lange     for (p = 0; p < numLeaves; ++p) {
8187bffde78SMichael Lange       localPointsNew[p] = localPoints[p];
8197bffde78SMichael Lange       remotePointsNew[p].index = remotePoints[p].index;
8207bffde78SMichael Lange       remotePointsNew[p].rank  = remotePoints[p].rank;
8217bffde78SMichael Lange     }
822f3190fc2SToby Isaac     p = numLeaves;
823e8f14785SLisandro Dalcin     ierr = PetscHMapIGetKeys(claimshash, &p, localPointsNew);CHKERRQ(ierr);
824f3190fc2SToby Isaac     ierr = PetscSortInt(numLocalNew, &localPointsNew[numLeaves]);CHKERRQ(ierr);
8257bffde78SMichael Lange     for (p = numLeaves; p < numLeaves + numLocalNew; ++p) {
826e8f14785SLisandro Dalcin       ierr = PetscHMapIGet(claimshash, localPointsNew[p], &offset);CHKERRQ(ierr);
8277bffde78SMichael Lange       remotePointsNew[p] = claims[offset];
8287bffde78SMichael Lange     }
8297bffde78SMichael Lange     ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), &sfPointNew);CHKERRQ(ierr);
8307bffde78SMichael Lange     ierr = PetscSFSetGraph(sfPointNew, pEnd-pStart, numLeaves+numLocalNew, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr);
8317bffde78SMichael Lange     ierr = DMSetPointSF(dm, sfPointNew);CHKERRQ(ierr);
8327bffde78SMichael Lange     ierr = PetscSFDestroy(&sfPointNew);CHKERRQ(ierr);
833e8f14785SLisandro Dalcin     ierr = PetscHMapIDestroy(&claimshash);CHKERRQ(ierr);
8347bffde78SMichael Lange   }
835e8f14785SLisandro Dalcin   ierr = PetscHMapIDestroy(&leafhash);CHKERRQ(ierr);
836e8f14785SLisandro Dalcin   ierr = PetscHMapIJDestroy(&roothash);CHKERRQ(ierr);
8377bffde78SMichael Lange   ierr = PetscSectionDestroy(&candidateSection);CHKERRQ(ierr);
8387bffde78SMichael Lange   ierr = PetscSectionDestroy(&candidateSectionRemote);CHKERRQ(ierr);
8397bffde78SMichael Lange   ierr = PetscSectionDestroy(&claimSection);CHKERRQ(ierr);
8407bffde78SMichael Lange   ierr = PetscFree(candidates);CHKERRQ(ierr);
8417bffde78SMichael Lange   ierr = PetscFree(candidatesRemote);CHKERRQ(ierr);
8427bffde78SMichael Lange   ierr = PetscFree(claims);CHKERRQ(ierr);
84325afeb17SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr);
8447bffde78SMichael Lange   PetscFunctionReturn(0);
8457bffde78SMichael Lange }
8467bffde78SMichael Lange 
8470c795ddcSMatthew G. Knepley /*@C
84880330477SMatthew G. Knepley   DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc.
84980330477SMatthew G. Knepley 
85080330477SMatthew G. Knepley   Collective on DM
85180330477SMatthew G. Knepley 
852e56d480eSMatthew G. Knepley   Input Parameters:
853e56d480eSMatthew G. Knepley + dm - The DMPlex object with only cells and vertices
85410a67516SMatthew G. Knepley - dmInt - The interpolated DM
85580330477SMatthew G. Knepley 
85680330477SMatthew G. Knepley   Output Parameter:
8574e3744c5SMatthew G. Knepley . dmInt - The complete DMPlex object
85880330477SMatthew G. Knepley 
85980330477SMatthew G. Knepley   Level: intermediate
86080330477SMatthew G. Knepley 
86195452b02SPatrick Sanan   Notes:
86295452b02SPatrick Sanan     It does not copy over the coordinates.
86343eeeb2dSStefano Zampini 
86480330477SMatthew G. Knepley .keywords: mesh
86543eeeb2dSStefano Zampini .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates()
86680330477SMatthew G. Knepley @*/
8679f074e33SMatthew G Knepley PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt)
8689f074e33SMatthew G Knepley {
8699a5b13a2SMatthew G Knepley   DM             idm, odm = dm;
8707bffde78SMichael Lange   PetscSF        sfPoint;
8712e1b13c2SMatthew G. Knepley   PetscInt       depth, dim, d;
87210a67516SMatthew G. Knepley   const char    *name;
8739f074e33SMatthew G Knepley   PetscErrorCode ierr;
8749f074e33SMatthew G Knepley 
8759f074e33SMatthew G Knepley   PetscFunctionBegin;
87610a67516SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
87710a67516SMatthew G. Knepley   PetscValidPointer(dmInt, 2);
878a72f3261SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr);
8792e1b13c2SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
880c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
881b21b8912SMatthew G. Knepley   if ((depth == dim) || (dim <= 1)) {
88276b791aaSMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
88376b791aaSMatthew G. Knepley     idm  = dm;
884b21b8912SMatthew G. Knepley   } else {
8859a5b13a2SMatthew G Knepley     for (d = 1; d < dim; ++d) {
8869a5b13a2SMatthew G Knepley       /* Create interpolated mesh */
88710a67516SMatthew G. Knepley       ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr);
8889a5b13a2SMatthew G Knepley       ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr);
889c73cfb54SMatthew G. Knepley       ierr = DMSetDimension(idm, dim);CHKERRQ(ierr);
8907bffde78SMichael Lange       if (depth > 0) {
8917bffde78SMichael Lange         ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr);
8927bffde78SMichael Lange         ierr = DMGetPointSF(odm, &sfPoint);CHKERRQ(ierr);
8937bffde78SMichael Lange         ierr = DMPlexInterpolatePointSF(idm, sfPoint, depth);CHKERRQ(ierr);
8947bffde78SMichael Lange       }
8959a5b13a2SMatthew G Knepley       if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);}
8969a5b13a2SMatthew G Knepley       odm = idm;
8979f074e33SMatthew G Knepley     }
89810a67516SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) dm,  &name);CHKERRQ(ierr);
89910a67516SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) idm,  name);CHKERRQ(ierr);
90010a67516SMatthew G. Knepley     ierr = DMPlexCopyCoordinates(dm, idm);CHKERRQ(ierr);
90110a67516SMatthew G. Knepley     ierr = DMCopyLabels(dm, idm);CHKERRQ(ierr);
90284699f70SSatish Balay   }
90343eeeb2dSStefano Zampini   {
90443eeeb2dSStefano Zampini     PetscBool            isper;
90543eeeb2dSStefano Zampini     const PetscReal      *maxCell, *L;
90643eeeb2dSStefano Zampini     const DMBoundaryType *bd;
90743eeeb2dSStefano Zampini 
90843eeeb2dSStefano Zampini     ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr);
90943eeeb2dSStefano Zampini     ierr = DMSetPeriodicity(idm,isper,maxCell,L,bd);CHKERRQ(ierr);
91043eeeb2dSStefano Zampini   }
9119a5b13a2SMatthew G Knepley   *dmInt = idm;
912a72f3261SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr);
9139f074e33SMatthew G Knepley   PetscFunctionReturn(0);
9149f074e33SMatthew G Knepley }
91507b0a1fcSMatthew G Knepley 
91680330477SMatthew G. Knepley /*@
91780330477SMatthew G. Knepley   DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices
91880330477SMatthew G. Knepley 
91980330477SMatthew G. Knepley   Collective on DM
92080330477SMatthew G. Knepley 
92180330477SMatthew G. Knepley   Input Parameter:
92280330477SMatthew G. Knepley . dmA - The DMPlex object with initial coordinates
92380330477SMatthew G. Knepley 
92480330477SMatthew G. Knepley   Output Parameter:
92580330477SMatthew G. Knepley . dmB - The DMPlex object with copied coordinates
92680330477SMatthew G. Knepley 
92780330477SMatthew G. Knepley   Level: intermediate
92880330477SMatthew G. Knepley 
92980330477SMatthew G. Knepley   Note: This is typically used when adding pieces other than vertices to a mesh
93080330477SMatthew G. Knepley 
93180330477SMatthew G. Knepley .keywords: mesh
93265f90189SJed Brown .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMGetCoordinateSection()
93380330477SMatthew G. Knepley @*/
93407b0a1fcSMatthew G Knepley PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB)
93507b0a1fcSMatthew G Knepley {
93607b0a1fcSMatthew G Knepley   Vec            coordinatesA, coordinatesB;
93743eeeb2dSStefano Zampini   VecType        vtype;
93807b0a1fcSMatthew G Knepley   PetscSection   coordSectionA, coordSectionB;
93907b0a1fcSMatthew G Knepley   PetscScalar   *coordsA, *coordsB;
9400bedd6beSMatthew G. Knepley   PetscInt       spaceDim, Nf, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d;
94143eeeb2dSStefano Zampini   PetscInt       cStartA, cEndA, cStartB, cEndB, cS, cE;
94243eeeb2dSStefano Zampini   PetscBool      lc = PETSC_FALSE;
94307b0a1fcSMatthew G Knepley   PetscErrorCode ierr;
94407b0a1fcSMatthew G Knepley 
94507b0a1fcSMatthew G Knepley   PetscFunctionBegin;
94643eeeb2dSStefano Zampini   PetscValidHeaderSpecific(dmA, DM_CLASSID, 1);
94743eeeb2dSStefano Zampini   PetscValidHeaderSpecific(dmB, DM_CLASSID, 2);
94876b791aaSMatthew G. Knepley   if (dmA == dmB) PetscFunctionReturn(0);
94907b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr);
95007b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr);
95107b0a1fcSMatthew 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);
95243eeeb2dSStefano Zampini   ierr = DMPlexGetHeightStratum(dmA, 0, &cStartA, &cEndA);CHKERRQ(ierr);
95343eeeb2dSStefano Zampini   ierr = DMPlexGetHeightStratum(dmB, 0, &cStartB, &cEndB);CHKERRQ(ierr);
95469d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr);
95569d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr);
956972bc18aSToby Isaac   if (coordSectionA == coordSectionB) PetscFunctionReturn(0);
9570bedd6beSMatthew G. Knepley   ierr = PetscSectionGetNumFields(coordSectionA, &Nf);CHKERRQ(ierr);
9580bedd6beSMatthew G. Knepley   if (!Nf) PetscFunctionReturn(0);
9590bedd6beSMatthew G. Knepley   if (Nf > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of coordinate fields must be 1, not %D", Nf);
960df26b574SMatthew G. Knepley   if (!coordSectionB) {
961df26b574SMatthew G. Knepley     PetscInt dim;
962df26b574SMatthew G. Knepley 
963df26b574SMatthew G. Knepley     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) coordSectionA), &coordSectionB);CHKERRQ(ierr);
964df26b574SMatthew G. Knepley     ierr = DMGetCoordinateDim(dmA, &dim);CHKERRQ(ierr);
965df26b574SMatthew G. Knepley     ierr = DMSetCoordinateSection(dmB, dim, coordSectionB);CHKERRQ(ierr);
966df26b574SMatthew G. Knepley     ierr = PetscObjectDereference((PetscObject) coordSectionB);CHKERRQ(ierr);
967df26b574SMatthew G. Knepley   }
96807b0a1fcSMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr);
96907b0a1fcSMatthew G Knepley   ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr);
97007b0a1fcSMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr);
97143eeeb2dSStefano Zampini   ierr = PetscSectionGetChart(coordSectionA, &cS, &cE);CHKERRQ(ierr);
97243eeeb2dSStefano Zampini   if (cStartA <= cS && cS < cEndA) { /* localized coordinates */
97343eeeb2dSStefano Zampini     if ((cEndA-cStartA) != (cEndB-cStartB)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of cellls in first DM %d != %d in the second DM", cEndA-cStartA, cEndB-cStartB);
97443eeeb2dSStefano Zampini     cS = cS - cStartA + cStartB;
97543eeeb2dSStefano Zampini     cE = vEndB;
97643eeeb2dSStefano Zampini     lc = PETSC_TRUE;
97743eeeb2dSStefano Zampini   } else {
97843eeeb2dSStefano Zampini     cS = vStartB;
97943eeeb2dSStefano Zampini     cE = vEndB;
98043eeeb2dSStefano Zampini   }
98143eeeb2dSStefano Zampini   ierr = PetscSectionSetChart(coordSectionB, cS, cE);CHKERRQ(ierr);
98207b0a1fcSMatthew G Knepley   for (v = vStartB; v < vEndB; ++v) {
98307b0a1fcSMatthew G Knepley     ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr);
98407b0a1fcSMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr);
98507b0a1fcSMatthew G Knepley   }
98643eeeb2dSStefano Zampini   if (lc) { /* localized coordinates */
98743eeeb2dSStefano Zampini     PetscInt c;
98843eeeb2dSStefano Zampini 
98943eeeb2dSStefano Zampini     for (c = cS-cStartB; c < cEndB-cStartB; c++) {
99043eeeb2dSStefano Zampini       PetscInt dof;
99143eeeb2dSStefano Zampini 
99243eeeb2dSStefano Zampini       ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr);
99343eeeb2dSStefano Zampini       ierr = PetscSectionSetDof(coordSectionB, c + cStartB, dof);CHKERRQ(ierr);
99443eeeb2dSStefano Zampini       ierr = PetscSectionSetFieldDof(coordSectionB, c + cStartB, 0, dof);CHKERRQ(ierr);
99543eeeb2dSStefano Zampini     }
99643eeeb2dSStefano Zampini   }
99707b0a1fcSMatthew G Knepley   ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr);
99807b0a1fcSMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr);
99907b0a1fcSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr);
10008b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinatesB);CHKERRQ(ierr);
100107b0a1fcSMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr);
100207b0a1fcSMatthew G Knepley   ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr);
100343eeeb2dSStefano Zampini   ierr = VecGetBlockSize(coordinatesA, &d);CHKERRQ(ierr);
100443eeeb2dSStefano Zampini   ierr = VecSetBlockSize(coordinatesB, d);CHKERRQ(ierr);
100543eeeb2dSStefano Zampini   ierr = VecGetType(coordinatesA, &vtype);CHKERRQ(ierr);
100643eeeb2dSStefano Zampini   ierr = VecSetType(coordinatesB, vtype);CHKERRQ(ierr);
100707b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr);
100807b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr);
100907b0a1fcSMatthew G Knepley   for (v = 0; v < vEndB-vStartB; ++v) {
101043eeeb2dSStefano Zampini     PetscInt offA, offB;
101143eeeb2dSStefano Zampini 
101243eeeb2dSStefano Zampini     ierr = PetscSectionGetOffset(coordSectionA, v + vStartA, &offA);CHKERRQ(ierr);
101343eeeb2dSStefano Zampini     ierr = PetscSectionGetOffset(coordSectionB, v + vStartB, &offB);CHKERRQ(ierr);
101407b0a1fcSMatthew G Knepley     for (d = 0; d < spaceDim; ++d) {
101543eeeb2dSStefano Zampini       coordsB[offB+d] = coordsA[offA+d];
101643eeeb2dSStefano Zampini     }
101743eeeb2dSStefano Zampini   }
101843eeeb2dSStefano Zampini   if (lc) { /* localized coordinates */
101943eeeb2dSStefano Zampini     PetscInt c;
102043eeeb2dSStefano Zampini 
102143eeeb2dSStefano Zampini     for (c = cS-cStartB; c < cEndB-cStartB; c++) {
102243eeeb2dSStefano Zampini       PetscInt dof, offA, offB;
102343eeeb2dSStefano Zampini 
102443eeeb2dSStefano Zampini       ierr = PetscSectionGetOffset(coordSectionA, c + cStartA, &offA);CHKERRQ(ierr);
102543eeeb2dSStefano Zampini       ierr = PetscSectionGetOffset(coordSectionB, c + cStartB, &offB);CHKERRQ(ierr);
102643eeeb2dSStefano Zampini       ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr);
102743eeeb2dSStefano Zampini       ierr = PetscMemcpy(coordsB + offB,coordsA + offA,dof*sizeof(*coordsB));CHKERRQ(ierr);
102807b0a1fcSMatthew G Knepley     }
102907b0a1fcSMatthew G Knepley   }
103007b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr);
103107b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr);
103207b0a1fcSMatthew G Knepley   ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr);
103307b0a1fcSMatthew G Knepley   ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr);
103407b0a1fcSMatthew G Knepley   PetscFunctionReturn(0);
103507b0a1fcSMatthew G Knepley }
10365c386225SMatthew G. Knepley 
10374e3744c5SMatthew G. Knepley /*@
10384e3744c5SMatthew G. Knepley   DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh
10394e3744c5SMatthew G. Knepley 
10404e3744c5SMatthew G. Knepley   Collective on DM
10414e3744c5SMatthew G. Knepley 
10424e3744c5SMatthew G. Knepley   Input Parameter:
10434e3744c5SMatthew G. Knepley . dm - The complete DMPlex object
10444e3744c5SMatthew G. Knepley 
10454e3744c5SMatthew G. Knepley   Output Parameter:
10464e3744c5SMatthew G. Knepley . dmUnint - The DMPlex object with only cells and vertices
10474e3744c5SMatthew G. Knepley 
10484e3744c5SMatthew G. Knepley   Level: intermediate
10494e3744c5SMatthew G. Knepley 
105095452b02SPatrick Sanan   Notes:
105195452b02SPatrick Sanan     It does not copy over the coordinates.
105243eeeb2dSStefano Zampini 
10534e3744c5SMatthew G. Knepley .keywords: mesh
105443eeeb2dSStefano Zampini .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates()
10554e3744c5SMatthew G. Knepley @*/
10564e3744c5SMatthew G. Knepley PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint)
10574e3744c5SMatthew G. Knepley {
10584e3744c5SMatthew G. Knepley   DM             udm;
1059c9f63434SStefano Zampini   PetscInt       dim, vStart, vEnd, cStart, cEnd, cMax, c, maxConeSize = 0, *cone;
10604e3744c5SMatthew G. Knepley   PetscErrorCode ierr;
10614e3744c5SMatthew G. Knepley 
10624e3744c5SMatthew G. Knepley   PetscFunctionBegin;
106343eeeb2dSStefano Zampini   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
106443eeeb2dSStefano Zampini   PetscValidPointer(dmUnint, 2);
1065c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
10664e3744c5SMatthew G. Knepley   if (dim <= 1) {
10674e3744c5SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
1068595d4abbSMatthew G. Knepley     *dmUnint = dm;
1069595d4abbSMatthew G. Knepley     PetscFunctionReturn(0);
10704e3744c5SMatthew G. Knepley   }
1071595d4abbSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1072595d4abbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1073c9f63434SStefano Zampini   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
10744e3744c5SMatthew G. Knepley   ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr);
10754e3744c5SMatthew G. Knepley   ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr);
1076c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(udm, dim);CHKERRQ(ierr);
10774e3744c5SMatthew G. Knepley   ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr);
10784e3744c5SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1079595d4abbSMatthew G. Knepley     PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
10804e3744c5SMatthew G. Knepley 
10814e3744c5SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
10824e3744c5SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
10834e3744c5SMatthew G. Knepley       const PetscInt p = closure[cl];
10844e3744c5SMatthew G. Knepley 
10854e3744c5SMatthew G. Knepley       if ((p >= vStart) && (p < vEnd)) ++coneSize;
10864e3744c5SMatthew G. Knepley     }
10874e3744c5SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
10884e3744c5SMatthew G. Knepley     ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr);
1089595d4abbSMatthew G. Knepley     maxConeSize = PetscMax(maxConeSize, coneSize);
10904e3744c5SMatthew G. Knepley   }
10914e3744c5SMatthew G. Knepley   ierr = DMSetUp(udm);CHKERRQ(ierr);
1092785e854fSJed Brown   ierr = PetscMalloc1(maxConeSize, &cone);CHKERRQ(ierr);
10934e3744c5SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1094595d4abbSMatthew G. Knepley     PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
10954e3744c5SMatthew G. Knepley 
10964e3744c5SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
10974e3744c5SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
10984e3744c5SMatthew G. Knepley       const PetscInt p = closure[cl];
10994e3744c5SMatthew G. Knepley 
11004e3744c5SMatthew G. Knepley       if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p;
11014e3744c5SMatthew G. Knepley     }
11024e3744c5SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
11034e3744c5SMatthew G. Knepley     ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr);
11044e3744c5SMatthew G. Knepley   }
11054e3744c5SMatthew G. Knepley   ierr = PetscFree(cone);CHKERRQ(ierr);
1106c9f63434SStefano Zampini   ierr = DMPlexSetHybridBounds(udm, cMax, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
11074e3744c5SMatthew G. Knepley   ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr);
11084e3744c5SMatthew G. Knepley   ierr = DMPlexStratify(udm);CHKERRQ(ierr);
11095c7de58aSMatthew G. Knepley   /* Reduce SF */
11105c7de58aSMatthew G. Knepley   {
11115c7de58aSMatthew G. Knepley     PetscSF            sfPoint, sfPointUn;
11125c7de58aSMatthew G. Knepley     const PetscSFNode *remotePoints;
11135c7de58aSMatthew G. Knepley     const PetscInt    *localPoints;
11145c7de58aSMatthew G. Knepley     PetscSFNode       *remotePointsUn;
11155c7de58aSMatthew G. Knepley     PetscInt          *localPointsUn;
11165c7de58aSMatthew G. Knepley     PetscInt           vEnd, numRoots, numLeaves, l;
11175c7de58aSMatthew G. Knepley     PetscInt           numLeavesUn = 0, n = 0;
11185c7de58aSMatthew G. Knepley     PetscErrorCode     ierr;
11195c7de58aSMatthew G. Knepley 
11205c7de58aSMatthew G. Knepley     /* Get original SF information */
11215c7de58aSMatthew G. Knepley     ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
11225c7de58aSMatthew G. Knepley     ierr = DMGetPointSF(udm, &sfPointUn);CHKERRQ(ierr);
11235c7de58aSMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, 0, NULL, &vEnd);CHKERRQ(ierr);
11245c7de58aSMatthew G. Knepley     ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr);
11255c7de58aSMatthew G. Knepley     /* Allocate space for cells and vertices */
11265c7de58aSMatthew G. Knepley     for (l = 0; l < numLeaves; ++l) if (localPoints[l] < vEnd) numLeavesUn++;
11275c7de58aSMatthew G. Knepley     /* Fill in leaves */
11285c7de58aSMatthew G. Knepley     if (vEnd >= 0) {
11295c7de58aSMatthew G. Knepley       ierr = PetscMalloc1(numLeavesUn, &remotePointsUn);CHKERRQ(ierr);
11305c7de58aSMatthew G. Knepley       ierr = PetscMalloc1(numLeavesUn, &localPointsUn);CHKERRQ(ierr);
11315c7de58aSMatthew G. Knepley       for (l = 0; l < numLeaves; l++) {
11325c7de58aSMatthew G. Knepley         if (localPoints[l] < vEnd) {
11335c7de58aSMatthew G. Knepley           localPointsUn[n]        = localPoints[l];
11345c7de58aSMatthew G. Knepley           remotePointsUn[n].rank  = remotePoints[l].rank;
11355c7de58aSMatthew G. Knepley           remotePointsUn[n].index = remotePoints[l].index;
11365c7de58aSMatthew G. Knepley           ++n;
11375c7de58aSMatthew G. Knepley         }
11385c7de58aSMatthew G. Knepley       }
11395c7de58aSMatthew G. Knepley       if (n != numLeavesUn) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent number of leaves %d != %d", n, numLeavesUn);
11405c7de58aSMatthew G. Knepley       ierr = PetscSFSetGraph(sfPointUn, vEnd, numLeavesUn, localPointsUn, PETSC_OWN_POINTER, remotePointsUn, PETSC_OWN_POINTER);CHKERRQ(ierr);
11415c7de58aSMatthew G. Knepley     }
11425c7de58aSMatthew G. Knepley   }
114343eeeb2dSStefano Zampini   {
114443eeeb2dSStefano Zampini     PetscBool            isper;
114543eeeb2dSStefano Zampini     const PetscReal      *maxCell, *L;
114643eeeb2dSStefano Zampini     const DMBoundaryType *bd;
114743eeeb2dSStefano Zampini 
114843eeeb2dSStefano Zampini     ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr);
114943eeeb2dSStefano Zampini     ierr = DMSetPeriodicity(udm,isper,maxCell,L,bd);CHKERRQ(ierr);
115043eeeb2dSStefano Zampini   }
115143eeeb2dSStefano Zampini 
11524e3744c5SMatthew G. Knepley   *dmUnint = udm;
11534e3744c5SMatthew G. Knepley   PetscFunctionReturn(0);
11544e3744c5SMatthew G. Knepley }
1155