xref: /petsc/src/dm/impls/plex/plexinterpolate.c (revision fa01033e04656531b15c4b2fb2e67046475708de)
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:
7799836e0eSStefano 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:
10599836e0eSStefano 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:
15799836e0eSStefano 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:
16199836e0eSStefano 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 
16699836e0eSStefano Zampini /*
16799836e0eSStefano Zampini   DMPlexGetRawFacesHybrid_Internal - Gets groups of vertices that correspond to faces for the given cone using hybrid ordering (prisms)
16899836e0eSStefano Zampini */
16999836e0eSStefano Zampini static PetscErrorCode DMPlexGetRawFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[])
17099836e0eSStefano Zampini {
17199836e0eSStefano Zampini   PetscInt       *facesTmp;
17299836e0eSStefano Zampini   PetscInt        maxConeSize, maxSupportSize;
17399836e0eSStefano Zampini   PetscErrorCode  ierr;
17499836e0eSStefano Zampini 
17599836e0eSStefano Zampini   PetscFunctionBegin;
17699836e0eSStefano Zampini   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
17799836e0eSStefano Zampini   if (faces && coneSize) PetscValidIntPointer(cone,4);
17899836e0eSStefano Zampini   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
17999836e0eSStefano Zampini   if (faces) {ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), MPIU_INT, &facesTmp);CHKERRQ(ierr);}
18099836e0eSStefano Zampini   switch (dim) {
18199836e0eSStefano Zampini   case 1:
18299836e0eSStefano Zampini     switch (coneSize) {
18399836e0eSStefano Zampini     case 2:
18499836e0eSStefano Zampini       if (faces) {
18599836e0eSStefano Zampini         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
18699836e0eSStefano Zampini         *faces = facesTmp;
18799836e0eSStefano Zampini       }
18899836e0eSStefano Zampini       if (numFaces)     *numFaces = 2;
18999836e0eSStefano Zampini       if (numFacesNotH) *numFacesNotH = 2;
19099836e0eSStefano Zampini       if (faceSize)     *faceSize = 1;
19199836e0eSStefano Zampini       break;
19299836e0eSStefano Zampini     default:
19399836e0eSStefano Zampini       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
19499836e0eSStefano Zampini     }
19599836e0eSStefano Zampini     break;
19699836e0eSStefano Zampini   case 2:
19799836e0eSStefano Zampini     switch (coneSize) {
19899836e0eSStefano Zampini     case 4:
19999836e0eSStefano Zampini       if (faces) {
20099836e0eSStefano Zampini         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
20199836e0eSStefano Zampini         facesTmp[2] = cone[2]; facesTmp[3] = cone[3];
20299836e0eSStefano Zampini         facesTmp[4] = cone[0]; facesTmp[5] = cone[2];
20399836e0eSStefano Zampini         facesTmp[6] = cone[1]; facesTmp[7] = cone[3];
20499836e0eSStefano Zampini         *faces = facesTmp;
20599836e0eSStefano Zampini       }
20699836e0eSStefano Zampini       if (numFaces)     *numFaces = 4;
20799836e0eSStefano Zampini       if (numFacesNotH) *numFacesNotH = 2;
20899836e0eSStefano Zampini       if (faceSize)     *faceSize = 2;
20999836e0eSStefano Zampini       break;
21099836e0eSStefano Zampini     default:
21199836e0eSStefano Zampini       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
21299836e0eSStefano Zampini     }
21399836e0eSStefano Zampini     break;
21499836e0eSStefano Zampini   case 3:
21599836e0eSStefano Zampini     switch (coneSize) {
21699836e0eSStefano Zampini     case 6: /* triangular prism */
21799836e0eSStefano Zampini       if (faces) {
21899836e0eSStefano Zampini         facesTmp[0]  = cone[0]; facesTmp[1]  = cone[1]; facesTmp[2]  = cone[2]; facesTmp[3]  = -1;      /* Bottom */
21999836e0eSStefano Zampini         facesTmp[4]  = cone[3]; facesTmp[5]  = cone[4]; facesTmp[6]  = cone[5]; facesTmp[7]  = -1;      /* Top */
22099836e0eSStefano Zampini         facesTmp[8]  = cone[0]; facesTmp[9]  = cone[1]; facesTmp[10] = cone[3]; facesTmp[11] = cone[4]; /* Back left */
22199836e0eSStefano Zampini         facesTmp[12] = cone[1]; facesTmp[13] = cone[2]; facesTmp[14] = cone[4]; facesTmp[15] = cone[5]; /* Back right */
22299836e0eSStefano Zampini         facesTmp[16] = cone[2]; facesTmp[17] = cone[0]; facesTmp[18] = cone[5]; facesTmp[19] = cone[3]; /* Front */
22399836e0eSStefano Zampini         *faces = facesTmp;
22499836e0eSStefano Zampini       }
22599836e0eSStefano Zampini       if (numFaces)     *numFaces = 5;
22699836e0eSStefano Zampini       if (numFacesNotH) *numFacesNotH = 2;
22799836e0eSStefano Zampini       if (faceSize)     *faceSize = -4;
22899836e0eSStefano Zampini       break;
22999836e0eSStefano Zampini     default:
23099836e0eSStefano Zampini       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
23199836e0eSStefano Zampini     }
23299836e0eSStefano Zampini     break;
23399836e0eSStefano Zampini   default:
23499836e0eSStefano Zampini     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim);
23599836e0eSStefano Zampini   }
23699836e0eSStefano Zampini   PetscFunctionReturn(0);
23799836e0eSStefano Zampini }
23899836e0eSStefano Zampini 
23999836e0eSStefano Zampini static PetscErrorCode DMPlexRestoreRawFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[])
24099836e0eSStefano Zampini {
24199836e0eSStefano Zampini   PetscErrorCode  ierr;
24299836e0eSStefano Zampini 
24399836e0eSStefano Zampini   PetscFunctionBegin;
24499836e0eSStefano Zampini   if (faces) { ierr = DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faces);CHKERRQ(ierr); }
24599836e0eSStefano Zampini   PetscFunctionReturn(0);
24699836e0eSStefano Zampini }
24799836e0eSStefano Zampini 
24899836e0eSStefano Zampini static PetscErrorCode DMPlexGetFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[])
24999836e0eSStefano Zampini {
25099836e0eSStefano Zampini   const PetscInt *cone = NULL;
25199836e0eSStefano Zampini   PetscInt        coneSize;
25299836e0eSStefano Zampini   PetscErrorCode  ierr;
25399836e0eSStefano Zampini 
25499836e0eSStefano Zampini   PetscFunctionBegin;
25599836e0eSStefano Zampini   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
25699836e0eSStefano Zampini   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
25799836e0eSStefano Zampini   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
25899836e0eSStefano Zampini   ierr = DMPlexGetRawFacesHybrid_Internal(dm, dim, coneSize, cone, numFaces, numFacesNotH, faceSize, faces);CHKERRQ(ierr);
25999836e0eSStefano Zampini   PetscFunctionReturn(0);
26099836e0eSStefano Zampini }
26199836e0eSStefano 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;
26999836e0eSStefano Zampini   PetscInt       coneSizeH = 0, faceSizeAllH = 0, numCellFacesH = 0, faceH, pMax = -1, dim, outerloop;
27099836e0eSStefano 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   }
29199836e0eSStefano Zampini   cMax = fMax = eMax = vMax = PETSC_DETERMINE;
29299836e0eSStefano Zampini   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
29399836e0eSStefano Zampini   if (cellDim == dim) {
29499836e0eSStefano Zampini     ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
29599836e0eSStefano Zampini     pMax = cMax;
29699836e0eSStefano Zampini   } else if (cellDim == dim -1) {
29799836e0eSStefano Zampini     ierr = DMPlexGetHybridBounds(dm, &cMax, &fMax, NULL, NULL);CHKERRQ(ierr);
29899836e0eSStefano Zampini     pMax = fMax;
29999836e0eSStefano Zampini   }
30099836e0eSStefano Zampini   pMax = pMax < 0 ? pEnd[cellDepth] : pMax;
30199836e0eSStefano Zampini   if (pMax < pEnd[cellDepth]) {
30299836e0eSStefano Zampini     const PetscInt *cellFaces, *cone;
30399836e0eSStefano Zampini     PetscInt        numCellFacesT, faceSize, cf;
3049f074e33SMatthew G Knepley 
30599836e0eSStefano Zampini     ierr = DMPlexGetConeSize(dm, pMax, &coneSizeH);CHKERRQ(ierr);
30699836e0eSStefano Zampini     ierr = DMPlexGetCone(dm, pMax, &cone);CHKERRQ(ierr);
30799836e0eSStefano Zampini     ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSizeH, cone, &numCellFacesH, &numCellFacesT, &faceSize, &cellFaces);CHKERRQ(ierr);
30899836e0eSStefano Zampini     if (faceSize < 0) {
30999836e0eSStefano Zampini       PetscInt *sizes, minv, maxv;
31099836e0eSStefano Zampini 
31199836e0eSStefano Zampini       /* count vertices of hybrid and non-hybrid faces */
31299836e0eSStefano Zampini       ierr = PetscCalloc1(numCellFacesH, &sizes);CHKERRQ(ierr);
31399836e0eSStefano Zampini       for (cf = 0; cf < numCellFacesT; ++cf) { /* These are the non-hybrid faces */
31499836e0eSStefano Zampini         const PetscInt *cellFace = &cellFaces[-cf*faceSize];
31599836e0eSStefano Zampini         PetscInt       f;
31699836e0eSStefano Zampini 
31799836e0eSStefano Zampini         for (f = 0; f < -faceSize; ++f) sizes[cf] += (cellFace[f] >= 0 ? 1 : 0);
31899836e0eSStefano Zampini       }
31999836e0eSStefano Zampini       ierr = PetscSortInt(numCellFacesT, sizes);CHKERRQ(ierr);
32099836e0eSStefano Zampini       minv = sizes[0];
32199836e0eSStefano Zampini       maxv = sizes[PetscMax(numCellFacesT-1, 0)];
32299836e0eSStefano Zampini       if (minv != maxv) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Different number of vertices for non-hybrid face %D != %D", minv, maxv);
32399836e0eSStefano Zampini       faceSizeAll = minv;
32499836e0eSStefano Zampini       ierr = PetscMemzero(sizes, numCellFacesH*sizeof(PetscInt));CHKERRQ(ierr);
32599836e0eSStefano Zampini       for (cf = numCellFacesT; cf < numCellFacesH; ++cf) { /* These are the hybrid faces */
32699836e0eSStefano Zampini         const PetscInt *cellFace = &cellFaces[-cf*faceSize];
32799836e0eSStefano Zampini         PetscInt       f;
32899836e0eSStefano Zampini 
32999836e0eSStefano Zampini         for (f = 0; f < -faceSize; ++f) sizes[cf-numCellFacesT] += (cellFace[f] >= 0 ? 1 : 0);
33099836e0eSStefano Zampini       }
33199836e0eSStefano Zampini       ierr = PetscSortInt(numCellFacesH - numCellFacesT, sizes);CHKERRQ(ierr);
33299836e0eSStefano Zampini       minv = sizes[0];
33399836e0eSStefano Zampini       maxv = sizes[PetscMax(numCellFacesH - numCellFacesT-1, 0)];
33499836e0eSStefano Zampini       ierr = PetscFree(sizes);CHKERRQ(ierr);
33599836e0eSStefano Zampini       if (minv != maxv) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Different number of vertices for hybrid face %D != %D", minv, maxv);
33699836e0eSStefano Zampini       faceSizeAllH = minv;
33799836e0eSStefano Zampini     } else { /* the size of the faces in hybrid cells is the same */
33899836e0eSStefano Zampini       faceSizeAll = faceSizeAllH = faceSize;
33999836e0eSStefano Zampini     }
34099836e0eSStefano Zampini     ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, cone, &numCellFacesH, &numCellFacesT, &faceSize, &cellFaces);CHKERRQ(ierr);
34199836e0eSStefano Zampini   } else if (pEnd[cellDepth] > pStart[cellDepth]) {
34299836e0eSStefano Zampini     ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr);
34399836e0eSStefano Zampini   }
34499836e0eSStefano Zampini   if (faceSizeAll > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll);
34599836e0eSStefano Zampini 
34699836e0eSStefano Zampini   /* With hybrid grids, we first iterate on hybrid cells and start numbering the non-hybrid grids
34799836e0eSStefano Zampini      Then, faces for non-hybrid cells are numbered.
34899836e0eSStefano Zampini      This is to guarantee consistent orientations (all 0) of all the points in the cone of the hybrid cells */
34999836e0eSStefano Zampini   ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr);
35099836e0eSStefano Zampini   for (outerloop = 0, face = pStart[faceDepth]; outerloop < 2; outerloop++) {
35199836e0eSStefano Zampini     PetscInt start, end;
35299836e0eSStefano Zampini 
35399836e0eSStefano Zampini     start = outerloop == 0 ? pMax : pStart[cellDepth];
35499836e0eSStefano Zampini     end = outerloop == 0 ? pEnd[cellDepth] : pMax;
35599836e0eSStefano Zampini     for (c = start; c < end; ++c) {
35699836e0eSStefano Zampini       const PetscInt *cellFaces;
35799836e0eSStefano Zampini       PetscInt        numCellFaces, faceSize, faceSizeInc, cf;
35899836e0eSStefano Zampini 
35999836e0eSStefano 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);
36299836e0eSStefano Zampini       } else { /* Hybrid cell */
36399836e0eSStefano Zampini         const PetscInt *cone;
36499836e0eSStefano Zampini         PetscInt        numCellFacesN, coneSize;
36599836e0eSStefano Zampini 
36699836e0eSStefano Zampini         ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
36799836e0eSStefano Zampini         if (coneSize != coneSizeH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid coneSize %D != %D", coneSize, coneSizeH);
36899836e0eSStefano Zampini         ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
36999836e0eSStefano Zampini         ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
37099836e0eSStefano Zampini         if (numCellFaces != numCellFacesH) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected numCellFaces %D != %D for hybrid cell %D", numCellFaces, numCellFacesH, c);
37199836e0eSStefano Zampini         faceSize = PetscMax(faceSize, -faceSize);
37299836e0eSStefano Zampini         if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSize);
37399836e0eSStefano Zampini         numCellFaces = numCellFacesN; /* process only non-hybrid faces */
37499836e0eSStefano Zampini       }
37599836e0eSStefano Zampini       faceSizeInc = faceSize;
3769f074e33SMatthew G Knepley       for (cf = 0; cf < numCellFaces; ++cf) {
37799836e0eSStefano Zampini         const PetscInt   *cellFace = &cellFaces[cf*faceSizeInc];
37899836e0eSStefano Zampini         PetscInt          faceSizeH = faceSize;
3799a5b13a2SMatthew G Knepley         PetscHashIJKLKey  key;
380e8f14785SLisandro Dalcin         PetscHashIter     iter;
381e8f14785SLisandro Dalcin         PetscBool         missing;
3829f074e33SMatthew G Knepley 
38399836e0eSStefano Zampini         if (faceSizeInc == 2) {
3849a5b13a2SMatthew G Knepley           key.i = PetscMin(cellFace[0], cellFace[1]);
3859a5b13a2SMatthew G Knepley           key.j = PetscMax(cellFace[0], cellFace[1]);
38699836e0eSStefano Zampini           key.k = PETSC_MAX_INT;
38799836e0eSStefano Zampini           key.l = PETSC_MAX_INT;
3889a5b13a2SMatthew G Knepley         } else {
38999836e0eSStefano Zampini           key.i = cellFace[0];
39099836e0eSStefano Zampini           key.j = cellFace[1];
39199836e0eSStefano Zampini           key.k = cellFace[2];
39299836e0eSStefano 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         }
39599836e0eSStefano Zampini         /* this check is redundant for non-hybrid meshes */
39699836e0eSStefano 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       }
40099836e0eSStefano Zampini       if (c < pMax) {
401439ece16SMatthew G. Knepley         ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
40299836e0eSStefano Zampini       } else {
40399836e0eSStefano Zampini         ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, NULL, NULL, NULL, NULL, &cellFaces);CHKERRQ(ierr);
40499836e0eSStefano Zampini       }
40599836e0eSStefano Zampini     }
40699836e0eSStefano Zampini   }
40799836e0eSStefano Zampini   pEnd[faceDepth] = face;
40899836e0eSStefano Zampini 
40999836e0eSStefano Zampini   /* Second pass for hybrid meshes: number hybrid faces */
41099836e0eSStefano Zampini   for (c = pMax; c < pEnd[cellDepth]; ++c) {
41199836e0eSStefano Zampini     const PetscInt *cellFaces, *cone;
41299836e0eSStefano Zampini     PetscInt        numCellFaces, numCellFacesN, faceSize, cf, coneSize;
41399836e0eSStefano Zampini 
41499836e0eSStefano Zampini     ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
41599836e0eSStefano Zampini     ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
41699836e0eSStefano Zampini     ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
41799836e0eSStefano Zampini     if (numCellFaces != numCellFacesH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid numCellFaces %D != %D", numCellFaces, numCellFacesH);
41899836e0eSStefano Zampini     faceSize = PetscMax(faceSize, -faceSize);
41999836e0eSStefano Zampini     for (cf = numCellFacesN; cf < numCellFaces; ++cf) { /* These are the hybrid faces */
42099836e0eSStefano Zampini       const PetscInt   *cellFace = &cellFaces[cf*faceSize];
42199836e0eSStefano Zampini       PetscHashIJKLKey  key;
42299836e0eSStefano Zampini       PetscHashIter     iter;
42399836e0eSStefano Zampini       PetscBool         missing;
42499836e0eSStefano Zampini       PetscInt          faceSizeH = faceSize;
42599836e0eSStefano Zampini 
42699836e0eSStefano Zampini       if (faceSize == 2) {
42799836e0eSStefano Zampini         key.i = PetscMin(cellFace[0], cellFace[1]);
42899836e0eSStefano Zampini         key.j = PetscMax(cellFace[0], cellFace[1]);
42999836e0eSStefano Zampini         key.k = PETSC_MAX_INT;
43099836e0eSStefano Zampini         key.l = PETSC_MAX_INT;
43199836e0eSStefano Zampini       } else {
43299836e0eSStefano Zampini         key.i = cellFace[0];
43399836e0eSStefano Zampini         key.j = cellFace[1];
43499836e0eSStefano Zampini         key.k = cellFace[2];
43599836e0eSStefano Zampini         key.l = faceSize > 3 ? (cellFace[3] < 0 ? faceSizeH = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT;
43699836e0eSStefano Zampini         ierr  = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr);
43799836e0eSStefano Zampini       }
43899836e0eSStefano 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);
43999836e0eSStefano Zampini       ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr);
44099836e0eSStefano Zampini       if (missing) {ierr = PetscHashIJKLIterSet(faceTable, iter, face++);CHKERRQ(ierr);}
44199836e0eSStefano Zampini     }
44299836e0eSStefano Zampini     ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
44399836e0eSStefano Zampini   }
44499836e0eSStefano Zampini   faceH = face - pEnd[faceDepth];
44599836e0eSStefano Zampini   if (faceH) {
44699836e0eSStefano Zampini     if (fMax == PETSC_DETERMINE) fMax = pEnd[faceDepth];
44799836e0eSStefano Zampini     else if (eMax == PETSC_DETERMINE) eMax = pEnd[faceDepth];
44899836e0eSStefano 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 */
46399836e0eSStefano Zampini       for (p = pStart[d]; p < pEnd[d]-faceH; ++p) {
4649a5b13a2SMatthew G Knepley         ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr);
4659a5b13a2SMatthew G Knepley       }
46699836e0eSStefano Zampini       for (p = pEnd[d]-faceH; p < pEnd[d]; ++p) {
46799836e0eSStefano Zampini         ierr = DMPlexSetConeSize(idm, p, faceSizeAllH);CHKERRQ(ierr);
46899836e0eSStefano 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*/
47299836e0eSStefano Zampini         if (p < pMax) {
473a014044eSMatthew G Knepley           ierr = DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);CHKERRQ(ierr);
47499836e0eSStefano Zampini         } else {
47599836e0eSStefano Zampini           ierr = DMPlexGetFacesHybrid_Internal(dm, cellDim, p, &coneSize, NULL, NULL, NULL);CHKERRQ(ierr);
47699836e0eSStefano 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 */
48899836e0eSStefano 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   }
50199836e0eSStefano Zampini   for (outerloop = 0, face = pStart[faceDepth]; outerloop < 2; outerloop++) {
50299836e0eSStefano Zampini     PetscInt start, end;
5039f074e33SMatthew G Knepley 
50499836e0eSStefano Zampini     start = outerloop == 0 ? pMax : pStart[cellDepth];
50599836e0eSStefano Zampini     end = outerloop == 0 ? pEnd[cellDepth] : pMax;
50699836e0eSStefano Zampini     for (c = start; c < end; ++c) {
50799836e0eSStefano Zampini       const PetscInt *cellFaces;
50899836e0eSStefano Zampini       PetscInt        numCellFaces, faceSize, faceSizeInc, cf;
50999836e0eSStefano Zampini 
51099836e0eSStefano 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);
51399836e0eSStefano Zampini       } else {
51499836e0eSStefano Zampini         const PetscInt *cone;
51599836e0eSStefano Zampini         PetscInt        numCellFacesN, coneSize;
51699836e0eSStefano Zampini 
51799836e0eSStefano Zampini         ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
51899836e0eSStefano Zampini         if (coneSize != coneSizeH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid coneSize %D != %D", coneSize, coneSizeH);
51999836e0eSStefano Zampini         ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
52099836e0eSStefano Zampini         ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
52199836e0eSStefano Zampini         if (numCellFaces != numCellFacesH) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected numCellFaces %D != %D for hybrid cell %D", numCellFaces, numCellFacesH, c);
52299836e0eSStefano Zampini         faceSize = PetscMax(faceSize, -faceSize);
52399836e0eSStefano Zampini         if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSize);
52499836e0eSStefano Zampini         numCellFaces = numCellFacesN; /* process only non-hybrid faces */
52599836e0eSStefano Zampini       }
52699836e0eSStefano Zampini       faceSizeInc = faceSize;
5279f074e33SMatthew G Knepley       for (cf = 0; cf < numCellFaces; ++cf) {
52899836e0eSStefano Zampini         const PetscInt  *cellFace = &cellFaces[cf*faceSizeInc];
5299a5b13a2SMatthew G Knepley         PetscHashIJKLKey key;
530e8f14785SLisandro Dalcin         PetscHashIter    iter;
531e8f14785SLisandro Dalcin         PetscBool        missing;
5329f074e33SMatthew G Knepley 
53399836e0eSStefano Zampini         if (faceSizeInc == 2) {
5349a5b13a2SMatthew G Knepley           key.i = PetscMin(cellFace[0], cellFace[1]);
5359a5b13a2SMatthew G Knepley           key.j = PetscMax(cellFace[0], cellFace[1]);
53699836e0eSStefano Zampini           key.k = PETSC_MAX_INT;
53799836e0eSStefano 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];
54299836e0eSStefano Zampini           key.l = faceSizeInc > 3 ? (cellFace[3] < 0 ? faceSize = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT;
54399836e0eSStefano 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       }
57899836e0eSStefano Zampini       if (c < pMax) {
579439ece16SMatthew G. Knepley         ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
58099836e0eSStefano Zampini       } else {
58199836e0eSStefano Zampini         ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, NULL, NULL, NULL, NULL, &cellFaces);CHKERRQ(ierr);
5829f074e33SMatthew G Knepley       }
58399836e0eSStefano Zampini     }
58499836e0eSStefano Zampini   }
58599836e0eSStefano Zampini   /* Second pass for hybrid meshes: orient hybrid faces */
58699836e0eSStefano Zampini   for (c = pMax; c < pEnd[cellDepth]; ++c) {
58799836e0eSStefano Zampini     const PetscInt *cellFaces, *cone;
58899836e0eSStefano Zampini     PetscInt        numCellFaces, numCellFacesN, faceSize, cf, coneSize;
58999836e0eSStefano Zampini 
59099836e0eSStefano Zampini     ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
59199836e0eSStefano Zampini     ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
59299836e0eSStefano Zampini     ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
59399836e0eSStefano Zampini     if (numCellFaces != numCellFacesH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid numCellFaces %D != %D", numCellFaces, numCellFacesH);
59499836e0eSStefano Zampini     faceSize = PetscMax(faceSize, -faceSize);
59599836e0eSStefano Zampini     for (cf = numCellFacesN; cf < numCellFaces; ++cf) { /* These are the hybrid faces */
59699836e0eSStefano Zampini       const PetscInt   *cellFace = &cellFaces[cf*faceSize];
59799836e0eSStefano Zampini       PetscHashIJKLKey key;
59899836e0eSStefano Zampini       PetscHashIter    iter;
59999836e0eSStefano Zampini       PetscBool        missing;
60099836e0eSStefano Zampini       PetscInt         faceSizeH = faceSize;
60199836e0eSStefano Zampini 
60299836e0eSStefano Zampini       if (faceSize == 2) {
60399836e0eSStefano Zampini         key.i = PetscMin(cellFace[0], cellFace[1]);
60499836e0eSStefano Zampini         key.j = PetscMax(cellFace[0], cellFace[1]);
60599836e0eSStefano Zampini         key.k = PETSC_MAX_INT;
60699836e0eSStefano Zampini         key.l = PETSC_MAX_INT;
60799836e0eSStefano Zampini       } else {
60899836e0eSStefano Zampini         key.i = cellFace[0];
60999836e0eSStefano Zampini         key.j = cellFace[1];
61099836e0eSStefano Zampini         key.k = cellFace[2];
61199836e0eSStefano Zampini         key.l = faceSize > 3 ? (cellFace[3] < 0 ? faceSizeH = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT;
61299836e0eSStefano Zampini         ierr  = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr);
61399836e0eSStefano Zampini       }
61499836e0eSStefano 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);
61599836e0eSStefano Zampini       ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr);
61699836e0eSStefano Zampini       if (missing) {
61799836e0eSStefano Zampini         ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr);
61899836e0eSStefano Zampini         ierr = PetscHashIJKLIterSet(faceTable, iter, face);CHKERRQ(ierr);
61999836e0eSStefano Zampini         ierr = DMPlexInsertCone(idm, c, cf, face++);CHKERRQ(ierr);
62099836e0eSStefano Zampini       } else {
62199836e0eSStefano Zampini         const PetscInt *cone;
62299836e0eSStefano Zampini         PetscInt        coneSize, ornt, i, j, f;
62399836e0eSStefano Zampini 
62499836e0eSStefano Zampini         ierr = PetscHashIJKLIterGet(faceTable, iter, &f);CHKERRQ(ierr);
62599836e0eSStefano Zampini         ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr);
62699836e0eSStefano Zampini         /* Orient face: Do not allow reverse orientation at the first vertex */
62799836e0eSStefano Zampini         ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr);
62899836e0eSStefano Zampini         ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr);
62999836e0eSStefano 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);
63099836e0eSStefano Zampini         /* - First find the initial vertex */
63199836e0eSStefano Zampini         for (i = 0; i < faceSizeH; ++i) if (cellFace[0] == cone[i]) break;
63299836e0eSStefano Zampini         /* - Try forward comparison */
63399836e0eSStefano Zampini         for (j = 0; j < faceSizeH; ++j) if (cellFace[j] != cone[(i+j)%faceSizeH]) break;
63499836e0eSStefano Zampini         if (j == faceSizeH) {
63599836e0eSStefano Zampini           if ((faceSizeH == 2) && (i == 1)) ornt = -2;
63699836e0eSStefano Zampini           else                             ornt = i;
63799836e0eSStefano Zampini         } else {
63899836e0eSStefano Zampini           /* - Try backward comparison */
63999836e0eSStefano Zampini           for (j = 0; j < faceSizeH; ++j) if (cellFace[j] != cone[(i+faceSizeH-j)%faceSizeH]) break;
64099836e0eSStefano Zampini           if (j == faceSizeH) {
64199836e0eSStefano Zampini             if (i == 0) ornt = -faceSizeH;
64299836e0eSStefano Zampini             else        ornt = -i;
64399836e0eSStefano Zampini           } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation");
64499836e0eSStefano Zampini         }
64599836e0eSStefano Zampini         ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr);
64699836e0eSStefano Zampini       }
64799836e0eSStefano Zampini     }
64899836e0eSStefano Zampini     ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
64999836e0eSStefano Zampini   }
65099836e0eSStefano 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);
65499836e0eSStefano 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 
660*fa01033eSVaclav Hapla PetscErrorCode DMPlexOrientCell(DM dm, PetscInt p, PetscInt masterConeSize, const PetscInt masterCone[])
661*fa01033eSVaclav Hapla {
662*fa01033eSVaclav Hapla   PetscInt i, j, k, maxConeSize, coneSize, coneConeSize, supportSize, supportConeSize;
663*fa01033eSVaclav Hapla   PetscInt start0, start1, start;
664*fa01033eSVaclav Hapla   PetscBool reverse0, reverse1, reverse;
665*fa01033eSVaclav Hapla   PetscInt newornt;
666*fa01033eSVaclav Hapla   const PetscInt *cone, *support, *supportCone, *ornts;
667*fa01033eSVaclav Hapla   PetscInt *newcone, *newornts;
668*fa01033eSVaclav Hapla   PetscErrorCode ierr;
669*fa01033eSVaclav Hapla 
670*fa01033eSVaclav Hapla   PetscFunctionBegin;
671*fa01033eSVaclav Hapla   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
672*fa01033eSVaclav Hapla   if (!coneSize) PetscFunctionReturn(0); /* do nothing for points with no cone */
673*fa01033eSVaclav Hapla   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
674*fa01033eSVaclav Hapla   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
675*fa01033eSVaclav Hapla   ierr = DMPlexFixFaceOrientations_Orient_Private(coneSize, masterConeSize, masterCone, cone, &start1, &reverse1);CHKERRQ(ierr);
676*fa01033eSVaclav Hapla   if (!start1 && !reverse1) PetscFunctionReturn(0);
677*fa01033eSVaclav Hapla   /* permute p's cone and orientations */
678*fa01033eSVaclav Hapla   ierr = DMPlexGetConeOrientation(dm, p, &ornts);CHKERRQ(ierr);
679*fa01033eSVaclav Hapla   ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newcone);CHKERRQ(ierr);
680*fa01033eSVaclav Hapla   ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newornts);CHKERRQ(ierr);
681*fa01033eSVaclav Hapla   ierr = DMPlexFixFaceOrientations_Permute_Private(coneSize, cone, start1, reverse1, newcone);CHKERRQ(ierr);
682*fa01033eSVaclav Hapla   ierr = DMPlexFixFaceOrientations_Permute_Private(coneSize, ornts, start1, reverse1, newornts);CHKERRQ(ierr);
683*fa01033eSVaclav Hapla   /* if direction of p (face) is flipped, flip also p's cone points (edges) */
684*fa01033eSVaclav Hapla   if (reverse1) {
685*fa01033eSVaclav Hapla     for (i=0; i<coneSize; i++) {
686*fa01033eSVaclav Hapla       ierr = DMPlexGetConeSize(dm, cone[i], &coneConeSize);CHKERRQ(ierr);
687*fa01033eSVaclav Hapla       ierr = DMPlexFixFaceOrientations_Translate_Private(newornts[i], &start0, &reverse0);CHKERRQ(ierr);
688*fa01033eSVaclav Hapla       ierr = DMPlexFixFaceOrientations_Combine_Private(coneConeSize, start0, reverse0, 1, PETSC_FALSE, &start, &reverse);CHKERRQ(ierr);
689*fa01033eSVaclav Hapla       ierr = DMPlexFixFaceOrientations_TranslateBack_Private(coneConeSize, start, reverse, &newornts[i]);CHKERRQ(ierr);
690*fa01033eSVaclav Hapla     }
691*fa01033eSVaclav Hapla   }
692*fa01033eSVaclav Hapla   ierr = DMPlexSetConeOrientation(dm, p, newornts);CHKERRQ(ierr);
693*fa01033eSVaclav Hapla   /* fix oriention of p within cones of p's support points */
694*fa01033eSVaclav Hapla   ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr);
695*fa01033eSVaclav Hapla   ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr);
696*fa01033eSVaclav Hapla   for (j=0; j<supportSize; j++) {
697*fa01033eSVaclav Hapla     ierr = DMPlexGetCone(dm, support[j], &supportCone);CHKERRQ(ierr);
698*fa01033eSVaclav Hapla     ierr = DMPlexGetConeSize(dm, support[j], &supportConeSize);CHKERRQ(ierr);
699*fa01033eSVaclav Hapla     ierr = DMPlexGetConeOrientation(dm, support[j], &ornts);CHKERRQ(ierr);
700*fa01033eSVaclav Hapla     for (k=0; k<supportConeSize; k++) {
701*fa01033eSVaclav Hapla       if (supportCone[k] != p) continue;
702*fa01033eSVaclav Hapla       ierr = DMPlexFixFaceOrientations_Translate_Private(ornts[k], &start0, &reverse0);CHKERRQ(ierr);
703*fa01033eSVaclav Hapla       ierr = DMPlexFixFaceOrientations_Combine_Private(coneSize, start0, reverse0, start1, reverse1, &start, &reverse);CHKERRQ(ierr);
704*fa01033eSVaclav Hapla       ierr = DMPlexFixFaceOrientations_TranslateBack_Private(coneSize, start, reverse, &newornt);CHKERRQ(ierr);
705*fa01033eSVaclav Hapla       ierr = DMPlexInsertConeOrientation(dm, support[j], k, newornt);CHKERRQ(ierr);
706*fa01033eSVaclav Hapla     }
707*fa01033eSVaclav Hapla   }
708*fa01033eSVaclav Hapla   /* rewrite cone */
709*fa01033eSVaclav Hapla   ierr = DMPlexSetCone(dm, p, newcone);CHKERRQ(ierr);
710*fa01033eSVaclav Hapla   ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newcone);CHKERRQ(ierr);
711*fa01033eSVaclav Hapla   ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newornts);CHKERRQ(ierr);
712*fa01033eSVaclav Hapla   PetscFunctionReturn(0);
713*fa01033eSVaclav Hapla }
714*fa01033eSVaclav Hapla 
7152e745bdaSMatthew G. Knepley static PetscErrorCode DMPlexOrientPointSF_Internal(DM dm, PetscSF sf)
7162e745bdaSMatthew G. Knepley {
717cae7fe92SVaclav Hapla   PetscInt          (*roots)[2], (*leaves)[2];
7188a650c75SVaclav Hapla   PetscMPIInt       (*rootsRanks)[2], (*leavesRanks)[2];
7192e745bdaSMatthew G. Knepley   const PetscInt    *locals;
7202e745bdaSMatthew G. Knepley   const PetscSFNode *remotes;
7218a650c75SVaclav Hapla   PetscInt           nroots, nleaves, p, c;
722c744bbe1SVaclav Hapla   const PetscInt    *cone;
723c744bbe1SVaclav Hapla   PetscInt           coneSize, ind0, ind1;
7242e745bdaSMatthew G. Knepley   MPI_Comm           comm;
7252e745bdaSMatthew G. Knepley   PetscMPIInt        rank;
7262e745bdaSMatthew G. Knepley   PetscInt           debug = 0;
7272e745bdaSMatthew G. Knepley   PetscErrorCode     ierr;
7282e745bdaSMatthew G. Knepley 
7292e745bdaSMatthew G. Knepley   PetscFunctionBegin;
7302e745bdaSMatthew G. Knepley   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &locals, &remotes);CHKERRQ(ierr);
7313ede9f65SMatthew G. Knepley   if (nroots < 0) PetscFunctionReturn(0);
7328a650c75SVaclav Hapla   ierr = PetscMalloc4(nroots, &roots, nroots, &leaves, nroots, &rootsRanks, nroots, &leavesRanks);CHKERRQ(ierr);
7332e745bdaSMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
7342e745bdaSMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
7352e745bdaSMatthew G. Knepley   if (debug && rank == 0) {ierr = PetscSynchronizedPrintf(comm, "Roots\n");CHKERRQ(ierr);}
7369e24d8a0SVaclav Hapla   for (p = 0; p < nroots; ++p) {
7379e24d8a0SVaclav Hapla     ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
7389e24d8a0SVaclav Hapla     ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
7392e745bdaSMatthew G. Knepley     /* Translate all points to root numbering */
7408a650c75SVaclav Hapla     for (c = 0; c < 2; c++) {
7418a650c75SVaclav Hapla       if (coneSize > 1) {
7428a650c75SVaclav Hapla         ierr = PetscFindInt(cone[c], nleaves, locals, &ind0);CHKERRQ(ierr);
7438a650c75SVaclav Hapla         if (ind0 < 0) {
7448a650c75SVaclav Hapla           roots[p][c] = cone[c];
7458a650c75SVaclav Hapla           rootsRanks[p][c] = rank;
746c8148b97SVaclav Hapla         } else {
7478a650c75SVaclav Hapla           roots[p][c] = remotes[ind0].index;
7488a650c75SVaclav Hapla           rootsRanks[p][c] = remotes[ind0].rank;
7498a650c75SVaclav Hapla         }
7508a650c75SVaclav Hapla       } else {
7518a650c75SVaclav Hapla         roots[p][c] = -1;
7528a650c75SVaclav Hapla         rootsRanks[p][c] = -1;
7532e745bdaSMatthew G. Knepley       }
7542e745bdaSMatthew G. Knepley     }
7558a650c75SVaclav Hapla   }
7568a650c75SVaclav Hapla   if (debug) {
7578a650c75SVaclav Hapla     for (p = 0; p < nroots; ++p) {
7588a650c75SVaclav Hapla       ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
7598a650c75SVaclav Hapla       ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
7608a650c75SVaclav Hapla       if (coneSize > 1) {
7618a650c75SVaclav Hapla         ierr = PetscSynchronizedPrintf(comm, "[%d]  %D: cone=[%D %D] roots=[%D %D] rootsRanks=[%D %D]\n", rank, p, cone[0], cone[1], roots[p][0], roots[p][1], rootsRanks[p][0], rootsRanks[p][1]);CHKERRQ(ierr);
7628a650c75SVaclav Hapla       }
7638a650c75SVaclav Hapla     }
7648a650c75SVaclav Hapla     ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);
7658a650c75SVaclav Hapla   }
7669e24d8a0SVaclav Hapla   for (p = 0; p < nroots; ++p) {
7679e24d8a0SVaclav Hapla     leaves[p][0] = -2;
7689e24d8a0SVaclav Hapla     leaves[p][1] = -2;
769c8148b97SVaclav Hapla   }
7702e745bdaSMatthew G. Knepley   ierr = PetscSFBcastBegin(sf, MPIU_2INT, roots, leaves);CHKERRQ(ierr);
7718a650c75SVaclav Hapla   ierr = PetscSFBcastBegin(sf, MPI_2INT, rootsRanks, leavesRanks);CHKERRQ(ierr);
7722e745bdaSMatthew G. Knepley   ierr = PetscSFBcastEnd(sf, MPIU_2INT, roots, leaves);CHKERRQ(ierr);
7738a650c75SVaclav Hapla   ierr = PetscSFBcastEnd(sf, MPI_2INT, rootsRanks, leavesRanks);CHKERRQ(ierr);
774c8148b97SVaclav Hapla   if (debug) {ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);}
775c8148b97SVaclav Hapla   if (debug && rank == 0) {ierr = PetscSynchronizedPrintf(comm, "Referred leaves\n");CHKERRQ(ierr);}
7769e24d8a0SVaclav Hapla   for (p = 0; p < nroots; ++p) {
7779e24d8a0SVaclav Hapla     if (leaves[p][0] < 0) continue;
7789e24d8a0SVaclav Hapla     ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
7799e24d8a0SVaclav Hapla     ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
780ea211068SVaclav Hapla     if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d]  %D: cone=[%D %D] leaves=[%D %D] roots=[%D %D] leavesRanks=[%D %D] rootsRanks=[%D %D]\n", rank, p, cone[0], cone[1], leaves[p][0], leaves[p][1], roots[p][0], roots[p][1], leavesRanks[p][0], leavesRanks[p][1], rootsRanks[p][0], rootsRanks[p][1]);CHKERRQ(ierr);}
781ea211068SVaclav Hapla     if ((leaves[p][0] != roots[p][0]) || (leaves[p][1] != roots[p][1]) || (leavesRanks[p][0] != rootsRanks[p][0]) || (leavesRanks[p][0] != rootsRanks[p][0])) {
7822e745bdaSMatthew G. Knepley       /* TODO Generalize this to correct an arbitrary orientation */
7839e24d8a0SVaclav Hapla       ierr = DMPlexReverseCell(dm, p);CHKERRQ(ierr);
7849e24d8a0SVaclav Hapla       if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d]  reversed point %D\n", rank, p);CHKERRQ(ierr);}
78523aaf07bSVaclav Hapla     }
78623aaf07bSVaclav Hapla #if defined(PETSC_USE_DEBUG)
7876dc21c25SVaclav Hapla     ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
7886dc21c25SVaclav Hapla     ierr = PetscFindInt(cone[0], nleaves, locals, &ind0);CHKERRQ(ierr);
7896dc21c25SVaclav Hapla     ierr = PetscFindInt(cone[1], nleaves, locals, &ind1);CHKERRQ(ierr);
7909e24d8a0SVaclav Hapla     if (ind0 < 0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point SF contains %D but is missing cone point %D = %D", p, 0, cone[0]);
7919e24d8a0SVaclav Hapla     if (ind1 < 0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point SF contains %D but is missing cone point %D = %D", p, 1, cone[1]);
7929e24d8a0SVaclav Hapla     if (leaves[p][0]  != remotes[ind0].index)
7939e24d8a0SVaclav Hapla       SETERRQ9(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D Cone[%D] %D --> (%D, %D) != %D Cone[%D] SF Point (%D, %D)", p, 0, cone[0], remotes[ind0].rank, remotes[ind0].index, leaves[p][0], 0, remotes[p].rank, remotes[p].index);
7949e24d8a0SVaclav Hapla     if (leaves[p][1] != remotes[ind1].index)
7959e24d8a0SVaclav Hapla       SETERRQ9(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D Cone[%D] %D --> (%D, %D) != %D Cone[%D] SF Point (%D, %D)", p, 1, cone[1], remotes[ind1].rank, remotes[ind1].index, leaves[p][1], 1, remotes[p].rank, remotes[p].index);
79623aaf07bSVaclav Hapla #endif
7972e745bdaSMatthew G. Knepley   }
7982e745bdaSMatthew G. Knepley   if (debug) {ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);}
7998a650c75SVaclav Hapla   ierr = PetscFree4(roots, leaves, rootsRanks, leavesRanks);CHKERRQ(ierr);
8002e745bdaSMatthew G. Knepley   PetscFunctionReturn(0);
8012e745bdaSMatthew G. Knepley }
8022e745bdaSMatthew G. Knepley 
8037bffde78SMichael Lange /* This interpolates the PointSF in parallel following local interpolation */
8047bffde78SMichael Lange static PetscErrorCode DMPlexInterpolatePointSF(DM dm, PetscSF pointSF, PetscInt depth)
8057bffde78SMichael Lange {
8069852e123SBarry Smith   PetscMPIInt        size, rank;
8077bffde78SMichael Lange   PetscInt           p, c, d, dof, offset;
808ffd6914dSSatish Balay   PetscInt           numLeaves, numRoots, candidatesSize, candidatesRemoteSize;
8097bffde78SMichael Lange   const PetscInt    *localPoints;
8107bffde78SMichael Lange   const PetscSFNode *remotePoints;
8117bffde78SMichael Lange   PetscSFNode       *candidates, *candidatesRemote, *claims;
8127bffde78SMichael Lange   PetscSection       candidateSection, candidateSectionRemote, claimSection;
813e8f14785SLisandro Dalcin   PetscHMapI         leafhash;
814e8f14785SLisandro Dalcin   PetscHMapIJ        roothash;
8157bffde78SMichael Lange   PetscHashIJKey     key;
8167bffde78SMichael Lange   PetscErrorCode     ierr;
8177bffde78SMichael Lange 
8187bffde78SMichael Lange   PetscFunctionBegin;
8199852e123SBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr);
8207bffde78SMichael Lange   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
8217bffde78SMichael Lange   ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr);
8229852e123SBarry Smith   if (size < 2 || numRoots < 0) PetscFunctionReturn(0);
82325afeb17SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr);
8247bffde78SMichael Lange   /* Build hashes of points in the SF for efficient lookup */
825e8f14785SLisandro Dalcin   ierr = PetscHMapICreate(&leafhash);CHKERRQ(ierr);
826e8f14785SLisandro Dalcin   ierr = PetscHMapIJCreate(&roothash);CHKERRQ(ierr);
8277bffde78SMichael Lange   for (p = 0; p < numLeaves; ++p) {
828e8f14785SLisandro Dalcin     ierr = PetscHMapISet(leafhash, localPoints[p], p);CHKERRQ(ierr);
829e8f14785SLisandro Dalcin     key.i = remotePoints[p].index;
830e8f14785SLisandro Dalcin     key.j = remotePoints[p].rank;
831e8f14785SLisandro Dalcin     ierr = PetscHMapIJSet(roothash, key, p);CHKERRQ(ierr);
8327bffde78SMichael Lange   }
8337bffde78SMichael Lange   /* Build a section / SFNode array of candidate points in the single-level adjacency of leaves,
8347bffde78SMichael Lange      where each candidate is defined by the root entry for the other vertex that defines the edge. */
8357bffde78SMichael Lange   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSection);CHKERRQ(ierr);
8367bffde78SMichael Lange   ierr = PetscSectionSetChart(candidateSection, 0, numRoots);CHKERRQ(ierr);
8377bffde78SMichael Lange   {
838ffd6914dSSatish Balay     PetscInt leaf, root, idx, a, *adj = NULL;
8397bffde78SMichael Lange     for (p = 0; p < numLeaves; ++p) {
8407bffde78SMichael Lange       PetscInt adjSize = PETSC_DETERMINE;
8417bffde78SMichael Lange       ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr);
8427bffde78SMichael Lange       for (a = 0; a < adjSize; ++a) {
843e8f14785SLisandro Dalcin         ierr = PetscHMapIGet(leafhash, adj[a], &leaf);CHKERRQ(ierr);
8447bffde78SMichael Lange         if (leaf >= 0) {ierr = PetscSectionAddDof(candidateSection, localPoints[p], 1);CHKERRQ(ierr);}
8457bffde78SMichael Lange       }
8467bffde78SMichael Lange     }
8477bffde78SMichael Lange     ierr = PetscSectionSetUp(candidateSection);CHKERRQ(ierr);
8487bffde78SMichael Lange     ierr = PetscSectionGetStorageSize(candidateSection, &candidatesSize);CHKERRQ(ierr);
8497bffde78SMichael Lange     ierr = PetscMalloc1(candidatesSize, &candidates);CHKERRQ(ierr);
8507bffde78SMichael Lange     for (p = 0; p < numLeaves; ++p) {
8517bffde78SMichael Lange       PetscInt adjSize = PETSC_DETERMINE;
8527bffde78SMichael Lange       ierr = PetscSectionGetOffset(candidateSection, localPoints[p], &offset);CHKERRQ(ierr);
8537bffde78SMichael Lange       ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr);
8547bffde78SMichael Lange       for (idx = 0, a = 0; a < adjSize; ++a) {
855e8f14785SLisandro Dalcin         ierr = PetscHMapIGet(leafhash, adj[a], &root);CHKERRQ(ierr);
8567bffde78SMichael Lange         if (root >= 0) candidates[offset+idx++] = remotePoints[root];
8577bffde78SMichael Lange       }
8587bffde78SMichael Lange     }
8597bffde78SMichael Lange     ierr = PetscFree(adj);CHKERRQ(ierr);
8607bffde78SMichael Lange   }
8617bffde78SMichael Lange   /* Gather candidate section / array pair into the root partition via inverse(multi(pointSF)). */
8627bffde78SMichael Lange   {
8637bffde78SMichael Lange     PetscSF   sfMulti, sfInverse, sfCandidates;
8647bffde78SMichael Lange     PetscInt *remoteOffsets;
8657bffde78SMichael Lange     ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr);
8667bffde78SMichael Lange     ierr = PetscSFCreateInverseSF(sfMulti, &sfInverse);CHKERRQ(ierr);
8677bffde78SMichael Lange     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSectionRemote);CHKERRQ(ierr);
8687bffde78SMichael Lange     ierr = PetscSFDistributeSection(sfInverse, candidateSection, &remoteOffsets, candidateSectionRemote);CHKERRQ(ierr);
8697bffde78SMichael Lange     ierr = PetscSFCreateSectionSF(sfInverse, candidateSection, remoteOffsets, candidateSectionRemote, &sfCandidates);CHKERRQ(ierr);
8707bffde78SMichael Lange     ierr = PetscSectionGetStorageSize(candidateSectionRemote, &candidatesRemoteSize);CHKERRQ(ierr);
8717bffde78SMichael Lange     ierr = PetscMalloc1(candidatesRemoteSize, &candidatesRemote);CHKERRQ(ierr);
8727bffde78SMichael Lange     ierr = PetscSFBcastBegin(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr);
8737bffde78SMichael Lange     ierr = PetscSFBcastEnd(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr);
8747bffde78SMichael Lange     ierr = PetscSFDestroy(&sfInverse);CHKERRQ(ierr);
8757bffde78SMichael Lange     ierr = PetscSFDestroy(&sfCandidates);CHKERRQ(ierr);
8767bffde78SMichael Lange     ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
8777bffde78SMichael Lange   }
8787bffde78SMichael Lange   /* Walk local roots and check for each remote candidate whether we know all required points,
8797bffde78SMichael Lange      either from owning it or having a root entry in the point SF. If we do we place a claim
8807bffde78SMichael Lange      by replacing the vertex number with our edge ID. */
8817bffde78SMichael Lange   {
8827bffde78SMichael Lange     PetscInt        idx, root, joinSize, vertices[2];
8837bffde78SMichael Lange     const PetscInt *rootdegree, *join = NULL;
8847bffde78SMichael Lange     ierr = PetscSFComputeDegreeBegin(pointSF, &rootdegree);CHKERRQ(ierr);
8857bffde78SMichael Lange     ierr = PetscSFComputeDegreeEnd(pointSF, &rootdegree);CHKERRQ(ierr);
8867bffde78SMichael Lange     /* Loop remote edge connections and put in a claim if both vertices are known */
8877bffde78SMichael Lange     for (idx = 0, p = 0; p < numRoots; ++p) {
8887bffde78SMichael Lange       for (d = 0; d < rootdegree[p]; ++d) {
8897bffde78SMichael Lange         ierr = PetscSectionGetDof(candidateSectionRemote, idx, &dof);CHKERRQ(ierr);
8907bffde78SMichael Lange         ierr = PetscSectionGetOffset(candidateSectionRemote, idx, &offset);CHKERRQ(ierr);
8917bffde78SMichael Lange         for (c = 0; c < dof; ++c) {
8927bffde78SMichael Lange           /* We own both vertices, so we claim the edge by replacing vertex with edge */
8937bffde78SMichael Lange           if (candidatesRemote[offset+c].rank == rank) {
8947bffde78SMichael Lange             vertices[0] = p; vertices[1] = candidatesRemote[offset+c].index;
8957bffde78SMichael Lange             ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
8967bffde78SMichael Lange             if (joinSize == 1) candidatesRemote[offset+c].index = join[0];
8977bffde78SMichael Lange             ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
8987bffde78SMichael Lange             continue;
8997bffde78SMichael Lange           }
9007bffde78SMichael Lange           /* If we own one vertex and share a root with the other, we claim it */
901e8f14785SLisandro Dalcin           key.i = candidatesRemote[offset+c].index;
902e8f14785SLisandro Dalcin           key.j = candidatesRemote[offset+c].rank;
903e8f14785SLisandro Dalcin           ierr = PetscHMapIJGet(roothash, key, &root);CHKERRQ(ierr);
9047bffde78SMichael Lange           if (root >= 0) {
9057bffde78SMichael Lange             vertices[0] = p; vertices[1] = localPoints[root];
9067bffde78SMichael Lange             ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
9077bffde78SMichael Lange             if (joinSize == 1) {
9087bffde78SMichael Lange               candidatesRemote[offset+c].index = join[0];
9097bffde78SMichael Lange               candidatesRemote[offset+c].rank = rank;
9107bffde78SMichael Lange             }
9117bffde78SMichael Lange             ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
9127bffde78SMichael Lange           }
9137bffde78SMichael Lange         }
9147bffde78SMichael Lange         idx++;
9157bffde78SMichael Lange       }
9167bffde78SMichael Lange     }
9177bffde78SMichael Lange   }
9187bffde78SMichael Lange   /* Push claims back to receiver via the MultiSF and derive new pointSF mapping on receiver */
9197bffde78SMichael Lange   {
9207bffde78SMichael Lange     PetscSF         sfMulti, sfClaims, sfPointNew;
921e8f14785SLisandro Dalcin     PetscHMapI      claimshash;
9227bffde78SMichael Lange     PetscInt        size, pStart, pEnd, root, joinSize, numLocalNew;
9237bffde78SMichael Lange     PetscInt       *remoteOffsets, *localPointsNew, vertices[2];
924ffd6914dSSatish Balay     const PetscInt *join = NULL;
9257bffde78SMichael Lange     PetscSFNode    *remotePointsNew;
9267bffde78SMichael Lange     ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr);
9277bffde78SMichael Lange     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &claimSection);CHKERRQ(ierr);
9287bffde78SMichael Lange     ierr = PetscSFDistributeSection(sfMulti, candidateSectionRemote, &remoteOffsets, claimSection);CHKERRQ(ierr);
9297bffde78SMichael Lange     ierr = PetscSFCreateSectionSF(sfMulti, candidateSectionRemote, remoteOffsets, claimSection, &sfClaims);CHKERRQ(ierr);
9307bffde78SMichael Lange     ierr = PetscSectionGetStorageSize(claimSection, &size);CHKERRQ(ierr);
9317bffde78SMichael Lange     ierr = PetscMalloc1(size, &claims);CHKERRQ(ierr);
9327bffde78SMichael Lange     ierr = PetscSFBcastBegin(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr);
9337bffde78SMichael Lange     ierr = PetscSFBcastEnd(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr);
9347bffde78SMichael Lange     ierr = PetscSFDestroy(&sfClaims);CHKERRQ(ierr);
9357bffde78SMichael Lange     ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
9367bffde78SMichael Lange     /* Walk the original section of local supports and add an SF entry for each updated item */
937e8f14785SLisandro Dalcin     ierr = PetscHMapICreate(&claimshash);CHKERRQ(ierr);
9387bffde78SMichael Lange     for (p = 0; p < numRoots; ++p) {
9397bffde78SMichael Lange       ierr = PetscSectionGetDof(candidateSection, p, &dof);CHKERRQ(ierr);
9407bffde78SMichael Lange       ierr = PetscSectionGetOffset(candidateSection, p, &offset);CHKERRQ(ierr);
9417bffde78SMichael Lange       for (d = 0; d < dof; ++d) {
9427bffde78SMichael Lange         if (candidates[offset+d].index != claims[offset+d].index) {
943e8f14785SLisandro Dalcin           key.i = candidates[offset+d].index;
944e8f14785SLisandro Dalcin           key.j = candidates[offset+d].rank;
945e8f14785SLisandro Dalcin           ierr = PetscHMapIJGet(roothash, key, &root);CHKERRQ(ierr);
9467bffde78SMichael Lange           if (root >= 0) {
9477bffde78SMichael Lange             vertices[0] = p; vertices[1] = localPoints[root];
9487bffde78SMichael Lange             ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
949e8f14785SLisandro Dalcin             if (joinSize == 1) {ierr = PetscHMapISet(claimshash, join[0], offset+d);CHKERRQ(ierr);}
9507bffde78SMichael Lange             ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr);
9517bffde78SMichael Lange           }
9527bffde78SMichael Lange         }
9537bffde78SMichael Lange       }
9547bffde78SMichael Lange     }
9557bffde78SMichael Lange     /* Create new pointSF from hashed claims */
956e8f14785SLisandro Dalcin     ierr = PetscHMapIGetSize(claimshash, &numLocalNew);CHKERRQ(ierr);
9577bffde78SMichael Lange     ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
9587bffde78SMichael Lange     ierr = PetscMalloc1(numLeaves + numLocalNew, &localPointsNew);CHKERRQ(ierr);
9597bffde78SMichael Lange     ierr = PetscMalloc1(numLeaves + numLocalNew, &remotePointsNew);CHKERRQ(ierr);
9607bffde78SMichael Lange     for (p = 0; p < numLeaves; ++p) {
9617bffde78SMichael Lange       localPointsNew[p] = localPoints[p];
9627bffde78SMichael Lange       remotePointsNew[p].index = remotePoints[p].index;
9637bffde78SMichael Lange       remotePointsNew[p].rank  = remotePoints[p].rank;
9647bffde78SMichael Lange     }
965f3190fc2SToby Isaac     p = numLeaves;
966e8f14785SLisandro Dalcin     ierr = PetscHMapIGetKeys(claimshash, &p, localPointsNew);CHKERRQ(ierr);
967f3190fc2SToby Isaac     ierr = PetscSortInt(numLocalNew, &localPointsNew[numLeaves]);CHKERRQ(ierr);
9687bffde78SMichael Lange     for (p = numLeaves; p < numLeaves + numLocalNew; ++p) {
969e8f14785SLisandro Dalcin       ierr = PetscHMapIGet(claimshash, localPointsNew[p], &offset);CHKERRQ(ierr);
9707bffde78SMichael Lange       remotePointsNew[p] = claims[offset];
9717bffde78SMichael Lange     }
9727bffde78SMichael Lange     ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), &sfPointNew);CHKERRQ(ierr);
9737bffde78SMichael Lange     ierr = PetscSFSetGraph(sfPointNew, pEnd-pStart, numLeaves+numLocalNew, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr);
9747bffde78SMichael Lange     ierr = DMSetPointSF(dm, sfPointNew);CHKERRQ(ierr);
9757bffde78SMichael Lange     ierr = PetscSFDestroy(&sfPointNew);CHKERRQ(ierr);
976e8f14785SLisandro Dalcin     ierr = PetscHMapIDestroy(&claimshash);CHKERRQ(ierr);
9777bffde78SMichael Lange   }
978e8f14785SLisandro Dalcin   ierr = PetscHMapIDestroy(&leafhash);CHKERRQ(ierr);
979e8f14785SLisandro Dalcin   ierr = PetscHMapIJDestroy(&roothash);CHKERRQ(ierr);
9807bffde78SMichael Lange   ierr = PetscSectionDestroy(&candidateSection);CHKERRQ(ierr);
9817bffde78SMichael Lange   ierr = PetscSectionDestroy(&candidateSectionRemote);CHKERRQ(ierr);
9827bffde78SMichael Lange   ierr = PetscSectionDestroy(&claimSection);CHKERRQ(ierr);
9837bffde78SMichael Lange   ierr = PetscFree(candidates);CHKERRQ(ierr);
9847bffde78SMichael Lange   ierr = PetscFree(candidatesRemote);CHKERRQ(ierr);
9857bffde78SMichael Lange   ierr = PetscFree(claims);CHKERRQ(ierr);
98625afeb17SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr);
9877bffde78SMichael Lange   PetscFunctionReturn(0);
9887bffde78SMichael Lange }
9897bffde78SMichael Lange 
9900c795ddcSMatthew G. Knepley /*@C
99180330477SMatthew G. Knepley   DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc.
99280330477SMatthew G. Knepley 
99380330477SMatthew G. Knepley   Collective on DM
99480330477SMatthew G. Knepley 
995e56d480eSMatthew G. Knepley   Input Parameters:
996e56d480eSMatthew G. Knepley + dm - The DMPlex object with only cells and vertices
99710a67516SMatthew G. Knepley - dmInt - The interpolated DM
99880330477SMatthew G. Knepley 
99980330477SMatthew G. Knepley   Output Parameter:
10004e3744c5SMatthew G. Knepley . dmInt - The complete DMPlex object
100180330477SMatthew G. Knepley 
100280330477SMatthew G. Knepley   Level: intermediate
100380330477SMatthew G. Knepley 
100495452b02SPatrick Sanan   Notes:
100595452b02SPatrick Sanan     It does not copy over the coordinates.
100643eeeb2dSStefano Zampini 
100780330477SMatthew G. Knepley .keywords: mesh
100843eeeb2dSStefano Zampini .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates()
100980330477SMatthew G. Knepley @*/
10109f074e33SMatthew G Knepley PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt)
10119f074e33SMatthew G Knepley {
10129a5b13a2SMatthew G Knepley   DM             idm, odm = dm;
10137bffde78SMichael Lange   PetscSF        sfPoint;
10142e1b13c2SMatthew G. Knepley   PetscInt       depth, dim, d;
101510a67516SMatthew G. Knepley   const char    *name;
10169f074e33SMatthew G Knepley   PetscErrorCode ierr;
10179f074e33SMatthew G Knepley 
10189f074e33SMatthew G Knepley   PetscFunctionBegin;
101910a67516SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
102010a67516SMatthew G. Knepley   PetscValidPointer(dmInt, 2);
1021a72f3261SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr);
10222e1b13c2SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1023c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1024b21b8912SMatthew G. Knepley   if ((depth == dim) || (dim <= 1)) {
102576b791aaSMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
102676b791aaSMatthew G. Knepley     idm  = dm;
1027b21b8912SMatthew G. Knepley   } else {
10289a5b13a2SMatthew G Knepley     for (d = 1; d < dim; ++d) {
10299a5b13a2SMatthew G Knepley       /* Create interpolated mesh */
103010a67516SMatthew G. Knepley       ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr);
10319a5b13a2SMatthew G Knepley       ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr);
1032c73cfb54SMatthew G. Knepley       ierr = DMSetDimension(idm, dim);CHKERRQ(ierr);
10337bffde78SMichael Lange       if (depth > 0) {
10347bffde78SMichael Lange         ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr);
10357bffde78SMichael Lange         ierr = DMGetPointSF(odm, &sfPoint);CHKERRQ(ierr);
10367bffde78SMichael Lange         ierr = DMPlexInterpolatePointSF(idm, sfPoint, depth);CHKERRQ(ierr);
10377bffde78SMichael Lange       }
10389a5b13a2SMatthew G Knepley       if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);}
10399a5b13a2SMatthew G Knepley       odm = idm;
10409f074e33SMatthew G Knepley     }
104110a67516SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) dm,  &name);CHKERRQ(ierr);
104210a67516SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) idm,  name);CHKERRQ(ierr);
104310a67516SMatthew G. Knepley     ierr = DMPlexCopyCoordinates(dm, idm);CHKERRQ(ierr);
104410a67516SMatthew G. Knepley     ierr = DMCopyLabels(dm, idm);CHKERRQ(ierr);
1045419e787bSVaclav Hapla     {
1046419e787bSVaclav Hapla       /* TODO temporary */
1047419e787bSVaclav Hapla       PetscBool flg=PETSC_FALSE;
1048419e787bSVaclav Hapla       ierr = PetscOptionsGetBool(NULL,NULL, "-dm_plex_check_point_sf", &flg, NULL);CHKERRQ(ierr);
1049419e787bSVaclav Hapla       if (flg) {ierr = DMPlexCheckPointSF(idm);CHKERRQ(ierr);}
1050419e787bSVaclav Hapla       ierr = DMViewFromOptions(idm, NULL, "-before_fix_dm_view");CHKERRQ(ierr);
1051419e787bSVaclav Hapla     }
10522e745bdaSMatthew G. Knepley     ierr = DMGetPointSF(idm, &sfPoint);CHKERRQ(ierr);
10532e745bdaSMatthew G. Knepley     ierr = DMPlexOrientPointSF_Internal(idm, sfPoint);CHKERRQ(ierr);
105484699f70SSatish Balay   }
105543eeeb2dSStefano Zampini   {
105643eeeb2dSStefano Zampini     PetscBool            isper;
105743eeeb2dSStefano Zampini     const PetscReal      *maxCell, *L;
105843eeeb2dSStefano Zampini     const DMBoundaryType *bd;
105943eeeb2dSStefano Zampini 
106043eeeb2dSStefano Zampini     ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr);
106143eeeb2dSStefano Zampini     ierr = DMSetPeriodicity(idm,isper,maxCell,L,bd);CHKERRQ(ierr);
106243eeeb2dSStefano Zampini   }
10639a5b13a2SMatthew G Knepley   *dmInt = idm;
1064a72f3261SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr);
10659f074e33SMatthew G Knepley   PetscFunctionReturn(0);
10669f074e33SMatthew G Knepley }
106707b0a1fcSMatthew G Knepley 
106880330477SMatthew G. Knepley /*@
106980330477SMatthew G. Knepley   DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices
107080330477SMatthew G. Knepley 
107180330477SMatthew G. Knepley   Collective on DM
107280330477SMatthew G. Knepley 
107380330477SMatthew G. Knepley   Input Parameter:
107480330477SMatthew G. Knepley . dmA - The DMPlex object with initial coordinates
107580330477SMatthew G. Knepley 
107680330477SMatthew G. Knepley   Output Parameter:
107780330477SMatthew G. Knepley . dmB - The DMPlex object with copied coordinates
107880330477SMatthew G. Knepley 
107980330477SMatthew G. Knepley   Level: intermediate
108080330477SMatthew G. Knepley 
108180330477SMatthew G. Knepley   Note: This is typically used when adding pieces other than vertices to a mesh
108280330477SMatthew G. Knepley 
108380330477SMatthew G. Knepley .keywords: mesh
108465f90189SJed Brown .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMGetCoordinateSection()
108580330477SMatthew G. Knepley @*/
108607b0a1fcSMatthew G Knepley PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB)
108707b0a1fcSMatthew G Knepley {
108807b0a1fcSMatthew G Knepley   Vec            coordinatesA, coordinatesB;
108943eeeb2dSStefano Zampini   VecType        vtype;
109007b0a1fcSMatthew G Knepley   PetscSection   coordSectionA, coordSectionB;
109107b0a1fcSMatthew G Knepley   PetscScalar   *coordsA, *coordsB;
10920bedd6beSMatthew G. Knepley   PetscInt       spaceDim, Nf, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d;
109343eeeb2dSStefano Zampini   PetscInt       cStartA, cEndA, cStartB, cEndB, cS, cE;
109443eeeb2dSStefano Zampini   PetscBool      lc = PETSC_FALSE;
109507b0a1fcSMatthew G Knepley   PetscErrorCode ierr;
109607b0a1fcSMatthew G Knepley 
109707b0a1fcSMatthew G Knepley   PetscFunctionBegin;
109843eeeb2dSStefano Zampini   PetscValidHeaderSpecific(dmA, DM_CLASSID, 1);
109943eeeb2dSStefano Zampini   PetscValidHeaderSpecific(dmB, DM_CLASSID, 2);
110076b791aaSMatthew G. Knepley   if (dmA == dmB) PetscFunctionReturn(0);
110107b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr);
110207b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr);
110307b0a1fcSMatthew 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);
110443eeeb2dSStefano Zampini   ierr = DMPlexGetHeightStratum(dmA, 0, &cStartA, &cEndA);CHKERRQ(ierr);
110543eeeb2dSStefano Zampini   ierr = DMPlexGetHeightStratum(dmB, 0, &cStartB, &cEndB);CHKERRQ(ierr);
110669d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr);
110769d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr);
1108972bc18aSToby Isaac   if (coordSectionA == coordSectionB) PetscFunctionReturn(0);
11090bedd6beSMatthew G. Knepley   ierr = PetscSectionGetNumFields(coordSectionA, &Nf);CHKERRQ(ierr);
11100bedd6beSMatthew G. Knepley   if (!Nf) PetscFunctionReturn(0);
11110bedd6beSMatthew G. Knepley   if (Nf > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of coordinate fields must be 1, not %D", Nf);
1112df26b574SMatthew G. Knepley   if (!coordSectionB) {
1113df26b574SMatthew G. Knepley     PetscInt dim;
1114df26b574SMatthew G. Knepley 
1115df26b574SMatthew G. Knepley     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) coordSectionA), &coordSectionB);CHKERRQ(ierr);
1116df26b574SMatthew G. Knepley     ierr = DMGetCoordinateDim(dmA, &dim);CHKERRQ(ierr);
1117df26b574SMatthew G. Knepley     ierr = DMSetCoordinateSection(dmB, dim, coordSectionB);CHKERRQ(ierr);
1118df26b574SMatthew G. Knepley     ierr = PetscObjectDereference((PetscObject) coordSectionB);CHKERRQ(ierr);
1119df26b574SMatthew G. Knepley   }
112007b0a1fcSMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr);
112107b0a1fcSMatthew G Knepley   ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr);
112207b0a1fcSMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr);
112343eeeb2dSStefano Zampini   ierr = PetscSectionGetChart(coordSectionA, &cS, &cE);CHKERRQ(ierr);
112443eeeb2dSStefano Zampini   if (cStartA <= cS && cS < cEndA) { /* localized coordinates */
112543eeeb2dSStefano 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);
112643eeeb2dSStefano Zampini     cS = cS - cStartA + cStartB;
112743eeeb2dSStefano Zampini     cE = vEndB;
112843eeeb2dSStefano Zampini     lc = PETSC_TRUE;
112943eeeb2dSStefano Zampini   } else {
113043eeeb2dSStefano Zampini     cS = vStartB;
113143eeeb2dSStefano Zampini     cE = vEndB;
113243eeeb2dSStefano Zampini   }
113343eeeb2dSStefano Zampini   ierr = PetscSectionSetChart(coordSectionB, cS, cE);CHKERRQ(ierr);
113407b0a1fcSMatthew G Knepley   for (v = vStartB; v < vEndB; ++v) {
113507b0a1fcSMatthew G Knepley     ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr);
113607b0a1fcSMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr);
113707b0a1fcSMatthew G Knepley   }
113843eeeb2dSStefano Zampini   if (lc) { /* localized coordinates */
113943eeeb2dSStefano Zampini     PetscInt c;
114043eeeb2dSStefano Zampini 
114143eeeb2dSStefano Zampini     for (c = cS-cStartB; c < cEndB-cStartB; c++) {
114243eeeb2dSStefano Zampini       PetscInt dof;
114343eeeb2dSStefano Zampini 
114443eeeb2dSStefano Zampini       ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr);
114543eeeb2dSStefano Zampini       ierr = PetscSectionSetDof(coordSectionB, c + cStartB, dof);CHKERRQ(ierr);
114643eeeb2dSStefano Zampini       ierr = PetscSectionSetFieldDof(coordSectionB, c + cStartB, 0, dof);CHKERRQ(ierr);
114743eeeb2dSStefano Zampini     }
114843eeeb2dSStefano Zampini   }
114907b0a1fcSMatthew G Knepley   ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr);
115007b0a1fcSMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr);
115107b0a1fcSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr);
11528b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinatesB);CHKERRQ(ierr);
115307b0a1fcSMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr);
115407b0a1fcSMatthew G Knepley   ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr);
115543eeeb2dSStefano Zampini   ierr = VecGetBlockSize(coordinatesA, &d);CHKERRQ(ierr);
115643eeeb2dSStefano Zampini   ierr = VecSetBlockSize(coordinatesB, d);CHKERRQ(ierr);
115743eeeb2dSStefano Zampini   ierr = VecGetType(coordinatesA, &vtype);CHKERRQ(ierr);
115843eeeb2dSStefano Zampini   ierr = VecSetType(coordinatesB, vtype);CHKERRQ(ierr);
115907b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr);
116007b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr);
116107b0a1fcSMatthew G Knepley   for (v = 0; v < vEndB-vStartB; ++v) {
116243eeeb2dSStefano Zampini     PetscInt offA, offB;
116343eeeb2dSStefano Zampini 
116443eeeb2dSStefano Zampini     ierr = PetscSectionGetOffset(coordSectionA, v + vStartA, &offA);CHKERRQ(ierr);
116543eeeb2dSStefano Zampini     ierr = PetscSectionGetOffset(coordSectionB, v + vStartB, &offB);CHKERRQ(ierr);
116607b0a1fcSMatthew G Knepley     for (d = 0; d < spaceDim; ++d) {
116743eeeb2dSStefano Zampini       coordsB[offB+d] = coordsA[offA+d];
116843eeeb2dSStefano Zampini     }
116943eeeb2dSStefano Zampini   }
117043eeeb2dSStefano Zampini   if (lc) { /* localized coordinates */
117143eeeb2dSStefano Zampini     PetscInt c;
117243eeeb2dSStefano Zampini 
117343eeeb2dSStefano Zampini     for (c = cS-cStartB; c < cEndB-cStartB; c++) {
117443eeeb2dSStefano Zampini       PetscInt dof, offA, offB;
117543eeeb2dSStefano Zampini 
117643eeeb2dSStefano Zampini       ierr = PetscSectionGetOffset(coordSectionA, c + cStartA, &offA);CHKERRQ(ierr);
117743eeeb2dSStefano Zampini       ierr = PetscSectionGetOffset(coordSectionB, c + cStartB, &offB);CHKERRQ(ierr);
117843eeeb2dSStefano Zampini       ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr);
117943eeeb2dSStefano Zampini       ierr = PetscMemcpy(coordsB + offB,coordsA + offA,dof*sizeof(*coordsB));CHKERRQ(ierr);
118007b0a1fcSMatthew G Knepley     }
118107b0a1fcSMatthew G Knepley   }
118207b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr);
118307b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr);
118407b0a1fcSMatthew G Knepley   ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr);
118507b0a1fcSMatthew G Knepley   ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr);
118607b0a1fcSMatthew G Knepley   PetscFunctionReturn(0);
118707b0a1fcSMatthew G Knepley }
11885c386225SMatthew G. Knepley 
11894e3744c5SMatthew G. Knepley /*@
11904e3744c5SMatthew G. Knepley   DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh
11914e3744c5SMatthew G. Knepley 
11924e3744c5SMatthew G. Knepley   Collective on DM
11934e3744c5SMatthew G. Knepley 
11944e3744c5SMatthew G. Knepley   Input Parameter:
11954e3744c5SMatthew G. Knepley . dm - The complete DMPlex object
11964e3744c5SMatthew G. Knepley 
11974e3744c5SMatthew G. Knepley   Output Parameter:
11984e3744c5SMatthew G. Knepley . dmUnint - The DMPlex object with only cells and vertices
11994e3744c5SMatthew G. Knepley 
12004e3744c5SMatthew G. Knepley   Level: intermediate
12014e3744c5SMatthew G. Knepley 
120295452b02SPatrick Sanan   Notes:
120395452b02SPatrick Sanan     It does not copy over the coordinates.
120443eeeb2dSStefano Zampini 
12054e3744c5SMatthew G. Knepley .keywords: mesh
120643eeeb2dSStefano Zampini .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates()
12074e3744c5SMatthew G. Knepley @*/
12084e3744c5SMatthew G. Knepley PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint)
12094e3744c5SMatthew G. Knepley {
12104e3744c5SMatthew G. Knepley   DM             udm;
1211c9f63434SStefano Zampini   PetscInt       dim, vStart, vEnd, cStart, cEnd, cMax, c, maxConeSize = 0, *cone;
12124e3744c5SMatthew G. Knepley   PetscErrorCode ierr;
12134e3744c5SMatthew G. Knepley 
12144e3744c5SMatthew G. Knepley   PetscFunctionBegin;
121543eeeb2dSStefano Zampini   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
121643eeeb2dSStefano Zampini   PetscValidPointer(dmUnint, 2);
1217c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
12184e3744c5SMatthew G. Knepley   if (dim <= 1) {
12194e3744c5SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
1220595d4abbSMatthew G. Knepley     *dmUnint = dm;
1221595d4abbSMatthew G. Knepley     PetscFunctionReturn(0);
12224e3744c5SMatthew G. Knepley   }
1223595d4abbSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1224595d4abbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1225c9f63434SStefano Zampini   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
12264e3744c5SMatthew G. Knepley   ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr);
12274e3744c5SMatthew G. Knepley   ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr);
1228c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(udm, dim);CHKERRQ(ierr);
12294e3744c5SMatthew G. Knepley   ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr);
12304e3744c5SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1231595d4abbSMatthew G. Knepley     PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
12324e3744c5SMatthew G. Knepley 
12334e3744c5SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
12344e3744c5SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
12354e3744c5SMatthew G. Knepley       const PetscInt p = closure[cl];
12364e3744c5SMatthew G. Knepley 
12374e3744c5SMatthew G. Knepley       if ((p >= vStart) && (p < vEnd)) ++coneSize;
12384e3744c5SMatthew G. Knepley     }
12394e3744c5SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
12404e3744c5SMatthew G. Knepley     ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr);
1241595d4abbSMatthew G. Knepley     maxConeSize = PetscMax(maxConeSize, coneSize);
12424e3744c5SMatthew G. Knepley   }
12434e3744c5SMatthew G. Knepley   ierr = DMSetUp(udm);CHKERRQ(ierr);
1244785e854fSJed Brown   ierr = PetscMalloc1(maxConeSize, &cone);CHKERRQ(ierr);
12454e3744c5SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1246595d4abbSMatthew G. Knepley     PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
12474e3744c5SMatthew G. Knepley 
12484e3744c5SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
12494e3744c5SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
12504e3744c5SMatthew G. Knepley       const PetscInt p = closure[cl];
12514e3744c5SMatthew G. Knepley 
12524e3744c5SMatthew G. Knepley       if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p;
12534e3744c5SMatthew G. Knepley     }
12544e3744c5SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
12554e3744c5SMatthew G. Knepley     ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr);
12564e3744c5SMatthew G. Knepley   }
12574e3744c5SMatthew G. Knepley   ierr = PetscFree(cone);CHKERRQ(ierr);
1258c9f63434SStefano Zampini   ierr = DMPlexSetHybridBounds(udm, cMax, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
12594e3744c5SMatthew G. Knepley   ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr);
12604e3744c5SMatthew G. Knepley   ierr = DMPlexStratify(udm);CHKERRQ(ierr);
12615c7de58aSMatthew G. Knepley   /* Reduce SF */
12625c7de58aSMatthew G. Knepley   {
12635c7de58aSMatthew G. Knepley     PetscSF            sfPoint, sfPointUn;
12645c7de58aSMatthew G. Knepley     const PetscSFNode *remotePoints;
12655c7de58aSMatthew G. Knepley     const PetscInt    *localPoints;
12665c7de58aSMatthew G. Knepley     PetscSFNode       *remotePointsUn;
12675c7de58aSMatthew G. Knepley     PetscInt          *localPointsUn;
12685c7de58aSMatthew G. Knepley     PetscInt           vEnd, numRoots, numLeaves, l;
12695c7de58aSMatthew G. Knepley     PetscInt           numLeavesUn = 0, n = 0;
12705c7de58aSMatthew G. Knepley     PetscErrorCode     ierr;
12715c7de58aSMatthew G. Knepley 
12725c7de58aSMatthew G. Knepley     /* Get original SF information */
12735c7de58aSMatthew G. Knepley     ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
12745c7de58aSMatthew G. Knepley     ierr = DMGetPointSF(udm, &sfPointUn);CHKERRQ(ierr);
12755c7de58aSMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, 0, NULL, &vEnd);CHKERRQ(ierr);
12765c7de58aSMatthew G. Knepley     ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr);
12775c7de58aSMatthew G. Knepley     /* Allocate space for cells and vertices */
12785c7de58aSMatthew G. Knepley     for (l = 0; l < numLeaves; ++l) if (localPoints[l] < vEnd) numLeavesUn++;
12795c7de58aSMatthew G. Knepley     /* Fill in leaves */
12805c7de58aSMatthew G. Knepley     if (vEnd >= 0) {
12815c7de58aSMatthew G. Knepley       ierr = PetscMalloc1(numLeavesUn, &remotePointsUn);CHKERRQ(ierr);
12825c7de58aSMatthew G. Knepley       ierr = PetscMalloc1(numLeavesUn, &localPointsUn);CHKERRQ(ierr);
12835c7de58aSMatthew G. Knepley       for (l = 0; l < numLeaves; l++) {
12845c7de58aSMatthew G. Knepley         if (localPoints[l] < vEnd) {
12855c7de58aSMatthew G. Knepley           localPointsUn[n]        = localPoints[l];
12865c7de58aSMatthew G. Knepley           remotePointsUn[n].rank  = remotePoints[l].rank;
12875c7de58aSMatthew G. Knepley           remotePointsUn[n].index = remotePoints[l].index;
12885c7de58aSMatthew G. Knepley           ++n;
12895c7de58aSMatthew G. Knepley         }
12905c7de58aSMatthew G. Knepley       }
12915c7de58aSMatthew G. Knepley       if (n != numLeavesUn) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent number of leaves %d != %d", n, numLeavesUn);
12925c7de58aSMatthew G. Knepley       ierr = PetscSFSetGraph(sfPointUn, vEnd, numLeavesUn, localPointsUn, PETSC_OWN_POINTER, remotePointsUn, PETSC_OWN_POINTER);CHKERRQ(ierr);
12935c7de58aSMatthew G. Knepley     }
12945c7de58aSMatthew G. Knepley   }
129543eeeb2dSStefano Zampini   {
129643eeeb2dSStefano Zampini     PetscBool            isper;
129743eeeb2dSStefano Zampini     const PetscReal      *maxCell, *L;
129843eeeb2dSStefano Zampini     const DMBoundaryType *bd;
129943eeeb2dSStefano Zampini 
130043eeeb2dSStefano Zampini     ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr);
130143eeeb2dSStefano Zampini     ierr = DMSetPeriodicity(udm,isper,maxCell,L,bd);CHKERRQ(ierr);
130243eeeb2dSStefano Zampini   }
130343eeeb2dSStefano Zampini 
13044e3744c5SMatthew G. Knepley   *dmUnint = udm;
13054e3744c5SMatthew G. Knepley   PetscFunctionReturn(0);
13064e3744c5SMatthew G. Knepley }
1307