xref: /petsc/src/dm/impls/plex/plexinterpolate.c (revision 66aa2a392c13f446a6cc6eae78cf84bbfd2ef47c)
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 
346b135d7daSStefano Zampini   /* With hybrid grids, we first iterate on hybrid cells and start numbering the non-hybrid faces
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 
660fa01033eSVaclav Hapla PetscErrorCode DMPlexOrientCell(DM dm, PetscInt p, PetscInt masterConeSize, const PetscInt masterCone[])
661fa01033eSVaclav Hapla {
66267aa0e3dSVaclav Hapla   PetscInt coneSize;
66367aa0e3dSVaclav Hapla   PetscInt start1=0;
66467aa0e3dSVaclav Hapla   PetscBool reverse1=PETSC_FALSE;
665a0d42dbcSVaclav Hapla   const PetscInt *cone=NULL;
666b8e38d09SVaclav Hapla   PetscErrorCode ierr;
667b8e38d09SVaclav Hapla 
668b8e38d09SVaclav Hapla   PetscFunctionBegin;
669b8e38d09SVaclav Hapla   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
670b8e38d09SVaclav Hapla   if (!coneSize) PetscFunctionReturn(0); /* do nothing for points with no cone */
671b8e38d09SVaclav Hapla   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
672b8e38d09SVaclav Hapla   ierr = DMPlexFixFaceOrientations_Orient_Private(coneSize, masterConeSize, masterCone, cone, &start1, &reverse1);CHKERRQ(ierr);
673b8e38d09SVaclav Hapla #if defined(PETSC_USE_DEBUG)
674b8e38d09SVaclav Hapla   if (PetscUnlikely(cone[start1] != masterCone[0])) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "The algorithm above is wrong as cone[%d] = %d != %d = masterCone[0]", start1, cone[start1], masterCone[0]);
675b8e38d09SVaclav Hapla #endif
676b8e38d09SVaclav Hapla   ierr = DMPlexOrientCell_Internal(dm, p, start1, reverse1);CHKERRQ(ierr);
677b8e38d09SVaclav Hapla #if defined(PETSC_USE_DEBUG)
67867aa0e3dSVaclav Hapla   {
67967aa0e3dSVaclav Hapla     PetscInt c;
680b8e38d09SVaclav Hapla     ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
681b8e38d09SVaclav Hapla     for (c = 0; c < 2; c++) {
682b8e38d09SVaclav Hapla       if (PetscUnlikely(cone[c] != masterCone[c])) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "The algorithm above is wrong as cone[%d] = %d != %d = masterCone[%d]", c, cone[c], masterCone[c], c);
683b8e38d09SVaclav Hapla     }
68467aa0e3dSVaclav Hapla   }
685b8e38d09SVaclav Hapla #endif
686b8e38d09SVaclav Hapla   PetscFunctionReturn(0);
687b8e38d09SVaclav Hapla }
688b8e38d09SVaclav Hapla 
6890fe46fbbSVaclav Hapla PetscErrorCode DMPlexOrientCell_Internal(DM dm, PetscInt p, PetscInt start1, PetscBool reverse1)
690b8e38d09SVaclav Hapla {
691fa01033eSVaclav Hapla   PetscInt i, j, k, maxConeSize, coneSize, coneConeSize, supportSize, supportConeSize;
692b8e38d09SVaclav Hapla   PetscInt start0, start;
693b8e38d09SVaclav Hapla   PetscBool reverse0, reverse;
694fa01033eSVaclav Hapla   PetscInt newornt;
695a0d42dbcSVaclav Hapla   const PetscInt *cone=NULL, *support=NULL, *supportCone=NULL, *ornts=NULL;
696a0d42dbcSVaclav Hapla   PetscInt *newcone=NULL, *newornts=NULL;
697fa01033eSVaclav Hapla   PetscErrorCode ierr;
698fa01033eSVaclav Hapla 
699fa01033eSVaclav Hapla   PetscFunctionBegin;
700b8e38d09SVaclav Hapla   if (!start1 && !reverse1) PetscFunctionReturn(0);
701fa01033eSVaclav Hapla   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
702fa01033eSVaclav Hapla   if (!coneSize) PetscFunctionReturn(0); /* do nothing for points with no cone */
703fa01033eSVaclav Hapla   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
704b8e38d09SVaclav Hapla   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
705fa01033eSVaclav Hapla   /* permute p's cone and orientations */
706fa01033eSVaclav Hapla   ierr = DMPlexGetConeOrientation(dm, p, &ornts);CHKERRQ(ierr);
707fa01033eSVaclav Hapla   ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newcone);CHKERRQ(ierr);
708fa01033eSVaclav Hapla   ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newornts);CHKERRQ(ierr);
709fa01033eSVaclav Hapla   ierr = DMPlexFixFaceOrientations_Permute_Private(coneSize, cone, start1, reverse1, newcone);CHKERRQ(ierr);
710fa01033eSVaclav Hapla   ierr = DMPlexFixFaceOrientations_Permute_Private(coneSize, ornts, start1, reverse1, newornts);CHKERRQ(ierr);
711fa01033eSVaclav Hapla   /* if direction of p (face) is flipped, flip also p's cone points (edges) */
712fa01033eSVaclav Hapla   if (reverse1) {
713fa01033eSVaclav Hapla     for (i=0; i<coneSize; i++) {
714fa01033eSVaclav Hapla       ierr = DMPlexGetConeSize(dm, cone[i], &coneConeSize);CHKERRQ(ierr);
715fa01033eSVaclav Hapla       ierr = DMPlexFixFaceOrientations_Translate_Private(newornts[i], &start0, &reverse0);CHKERRQ(ierr);
716fa01033eSVaclav Hapla       ierr = DMPlexFixFaceOrientations_Combine_Private(coneConeSize, start0, reverse0, 1, PETSC_FALSE, &start, &reverse);CHKERRQ(ierr);
717fa01033eSVaclav Hapla       ierr = DMPlexFixFaceOrientations_TranslateBack_Private(coneConeSize, start, reverse, &newornts[i]);CHKERRQ(ierr);
718fa01033eSVaclav Hapla     }
719fa01033eSVaclav Hapla   }
720fa01033eSVaclav Hapla   ierr = DMPlexSetConeOrientation(dm, p, newornts);CHKERRQ(ierr);
721fa01033eSVaclav Hapla   /* fix oriention of p within cones of p's support points */
722fa01033eSVaclav Hapla   ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr);
723fa01033eSVaclav Hapla   ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr);
724fa01033eSVaclav Hapla   for (j=0; j<supportSize; j++) {
725fa01033eSVaclav Hapla     ierr = DMPlexGetCone(dm, support[j], &supportCone);CHKERRQ(ierr);
726fa01033eSVaclav Hapla     ierr = DMPlexGetConeSize(dm, support[j], &supportConeSize);CHKERRQ(ierr);
727fa01033eSVaclav Hapla     ierr = DMPlexGetConeOrientation(dm, support[j], &ornts);CHKERRQ(ierr);
728fa01033eSVaclav Hapla     for (k=0; k<supportConeSize; k++) {
729fa01033eSVaclav Hapla       if (supportCone[k] != p) continue;
730fa01033eSVaclav Hapla       ierr = DMPlexFixFaceOrientations_Translate_Private(ornts[k], &start0, &reverse0);CHKERRQ(ierr);
731fa01033eSVaclav Hapla       ierr = DMPlexFixFaceOrientations_Combine_Private(coneSize, start0, reverse0, start1, reverse1, &start, &reverse);CHKERRQ(ierr);
732fa01033eSVaclav Hapla       ierr = DMPlexFixFaceOrientations_TranslateBack_Private(coneSize, start, reverse, &newornt);CHKERRQ(ierr);
733fa01033eSVaclav Hapla       ierr = DMPlexInsertConeOrientation(dm, support[j], k, newornt);CHKERRQ(ierr);
734fa01033eSVaclav Hapla     }
735fa01033eSVaclav Hapla   }
736fa01033eSVaclav Hapla   /* rewrite cone */
737fa01033eSVaclav Hapla   ierr = DMPlexSetCone(dm, p, newcone);CHKERRQ(ierr);
738fa01033eSVaclav Hapla   ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newcone);CHKERRQ(ierr);
739fa01033eSVaclav Hapla   ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newornts);CHKERRQ(ierr);
740fa01033eSVaclav Hapla   PetscFunctionReturn(0);
741fa01033eSVaclav Hapla }
742fa01033eSVaclav Hapla 
743f80536cbSVaclav Hapla static PetscErrorCode SortRmineRremoteByRemote_Private(PetscSF sf, PetscInt *rmine1[], PetscInt *rremote1[])
744f80536cbSVaclav Hapla {
745f80536cbSVaclav Hapla   PetscInt            nleaves;
746f80536cbSVaclav Hapla   PetscInt            nranks;
747a0d42dbcSVaclav Hapla   const PetscMPIInt  *ranks=NULL;
748a0d42dbcSVaclav Hapla   const PetscInt     *roffset=NULL, *rmine=NULL, *rremote=NULL;
749f80536cbSVaclav Hapla   PetscInt            n, o, r;
750f80536cbSVaclav Hapla   PetscErrorCode      ierr;
751f80536cbSVaclav Hapla 
752f80536cbSVaclav Hapla   PetscFunctionBegin;
753f80536cbSVaclav Hapla   ierr = PetscSFGetRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote);CHKERRQ(ierr);
754f80536cbSVaclav Hapla   nleaves = roffset[nranks];
755f80536cbSVaclav Hapla   ierr = PetscMalloc2(nleaves, rmine1, nleaves, rremote1);CHKERRQ(ierr);
756f80536cbSVaclav Hapla   for (r=0; r<nranks; r++) {
757f80536cbSVaclav Hapla     /* simultaneously sort rank-wise portions of rmine & rremote by values in rremote
758f80536cbSVaclav Hapla        - to unify order with the other side */
759f80536cbSVaclav Hapla     o = roffset[r];
760f80536cbSVaclav Hapla     n = roffset[r+1] - o;
761f80536cbSVaclav Hapla     ierr = PetscMemcpy(&(*rmine1)[o], &rmine[o], n*sizeof(PetscInt));CHKERRQ(ierr);
762f80536cbSVaclav Hapla     ierr = PetscMemcpy(&(*rremote1)[o], &rremote[o], n*sizeof(PetscInt));CHKERRQ(ierr);
763f80536cbSVaclav Hapla     ierr = PetscSortIntWithArray(n, &(*rremote1)[o], &(*rmine1)[o]);CHKERRQ(ierr);
764f80536cbSVaclav Hapla   }
765f80536cbSVaclav Hapla   PetscFunctionReturn(0);
766f80536cbSVaclav Hapla }
767f80536cbSVaclav Hapla 
768bba29ab8SVaclav Hapla PetscErrorCode DMPlexOrientInterface(DM dm)
7692e745bdaSMatthew G. Knepley {
770a0d42dbcSVaclav Hapla   PetscSF           sf=NULL;
771cae7fe92SVaclav Hapla   PetscInt          (*roots)[2], (*leaves)[2];
7728a650c75SVaclav Hapla   PetscMPIInt       (*rootsRanks)[2], (*leavesRanks)[2];
773a0d42dbcSVaclav Hapla   const PetscInt    *locals=NULL;
774a0d42dbcSVaclav Hapla   const PetscSFNode *remotes=NULL;
7758a650c75SVaclav Hapla   PetscInt           nroots, nleaves, p, c;
776f80536cbSVaclav Hapla   PetscInt           nranks, n, o, r;
777a0d42dbcSVaclav Hapla   const PetscMPIInt *ranks=NULL;
778a0d42dbcSVaclav Hapla   const PetscInt    *roffset=NULL;
779a0d42dbcSVaclav Hapla   PetscInt          *rmine1=NULL, *rremote1=NULL; /* rmine and rremote copies simultaneously sorted by rank and rremote */
780a0d42dbcSVaclav Hapla   const PetscInt    *cone=NULL;
781adeface4SVaclav Hapla   PetscInt           coneSize, ind0;
7822e745bdaSMatthew G. Knepley   MPI_Comm           comm;
7832e745bdaSMatthew G. Knepley   PetscMPIInt        rank;
7842e745bdaSMatthew G. Knepley   PetscInt           debug = 0;
7852e745bdaSMatthew G. Knepley   PetscErrorCode     ierr;
7862e745bdaSMatthew G. Knepley 
7872e745bdaSMatthew G. Knepley   PetscFunctionBegin;
788df6a9fadSVaclav Hapla   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
7892e745bdaSMatthew G. Knepley   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &locals, &remotes);CHKERRQ(ierr);
7903ede9f65SMatthew G. Knepley   if (nroots < 0) PetscFunctionReturn(0);
791f80536cbSVaclav Hapla   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
792f80536cbSVaclav Hapla   ierr = PetscSFGetRanks(sf, &nranks, &ranks, &roffset, NULL, NULL);CHKERRQ(ierr);
793dc21a0bfSVaclav Hapla #if defined(PETSC_USE_DEBUG)
794dc21a0bfSVaclav Hapla   ierr = DMViewFromOptions(dm, NULL, "-before_fix_dm_view");CHKERRQ(ierr);
795dc21a0bfSVaclav Hapla   ierr = DMPlexCheckPointSF(dm);CHKERRQ(ierr);
796dc21a0bfSVaclav Hapla #endif
797f80536cbSVaclav Hapla   ierr = SortRmineRremoteByRemote_Private(sf, &rmine1, &rremote1);CHKERRQ(ierr);
7988a650c75SVaclav Hapla   ierr = PetscMalloc4(nroots, &roots, nroots, &leaves, nroots, &rootsRanks, nroots, &leavesRanks);CHKERRQ(ierr);
7992e745bdaSMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
8002e745bdaSMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
8012e745bdaSMatthew G. Knepley   if (debug && rank == 0) {ierr = PetscSynchronizedPrintf(comm, "Roots\n");CHKERRQ(ierr);}
8029e24d8a0SVaclav Hapla   for (p = 0; p < nroots; ++p) {
8039e24d8a0SVaclav Hapla     ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
8049e24d8a0SVaclav Hapla     ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
8052e745bdaSMatthew G. Knepley     /* Translate all points to root numbering */
8068a650c75SVaclav Hapla     for (c = 0; c < 2; c++) {
8078a650c75SVaclav Hapla       if (coneSize > 1) {
8088a650c75SVaclav Hapla         ierr = PetscFindInt(cone[c], nleaves, locals, &ind0);CHKERRQ(ierr);
8098a650c75SVaclav Hapla         if (ind0 < 0) {
8108a650c75SVaclav Hapla           roots[p][c] = cone[c];
8118a650c75SVaclav Hapla           rootsRanks[p][c] = rank;
812c8148b97SVaclav Hapla         } else {
8138a650c75SVaclav Hapla           roots[p][c] = remotes[ind0].index;
8148a650c75SVaclav Hapla           rootsRanks[p][c] = remotes[ind0].rank;
8158a650c75SVaclav Hapla         }
8168a650c75SVaclav Hapla       } else {
8178a650c75SVaclav Hapla         roots[p][c] = -1;
8188a650c75SVaclav Hapla         rootsRanks[p][c] = -1;
8192e745bdaSMatthew G. Knepley       }
8202e745bdaSMatthew G. Knepley     }
8218a650c75SVaclav Hapla   }
8228a650c75SVaclav Hapla   if (debug) {
8238a650c75SVaclav Hapla     for (p = 0; p < nroots; ++p) {
8248a650c75SVaclav Hapla       ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
8258a650c75SVaclav Hapla       ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
8268a650c75SVaclav Hapla       if (coneSize > 1) {
8278a650c75SVaclav 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);
8288a650c75SVaclav Hapla       }
8298a650c75SVaclav Hapla     }
8308a650c75SVaclav Hapla     ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);
8318a650c75SVaclav Hapla   }
8329e24d8a0SVaclav Hapla   for (p = 0; p < nroots; ++p) {
8338ccb905dSVaclav Hapla     for (c = 0; c < 2; c++) {
8348ccb905dSVaclav Hapla       leaves[p][c] = -2;
8358ccb905dSVaclav Hapla       leavesRanks[p][c] = -2;
8368ccb905dSVaclav Hapla     }
837c8148b97SVaclav Hapla   }
8382e745bdaSMatthew G. Knepley   ierr = PetscSFBcastBegin(sf, MPIU_2INT, roots, leaves);CHKERRQ(ierr);
8398a650c75SVaclav Hapla   ierr = PetscSFBcastBegin(sf, MPI_2INT, rootsRanks, leavesRanks);CHKERRQ(ierr);
8402e745bdaSMatthew G. Knepley   ierr = PetscSFBcastEnd(sf, MPIU_2INT, roots, leaves);CHKERRQ(ierr);
8418a650c75SVaclav Hapla   ierr = PetscSFBcastEnd(sf, MPI_2INT, rootsRanks, leavesRanks);CHKERRQ(ierr);
842c8148b97SVaclav Hapla   if (debug) {ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);}
843c8148b97SVaclav Hapla   if (debug && rank == 0) {ierr = PetscSynchronizedPrintf(comm, "Referred leaves\n");CHKERRQ(ierr);}
8449e24d8a0SVaclav Hapla   for (p = 0; p < nroots; ++p) {
8459e24d8a0SVaclav Hapla     if (leaves[p][0] < 0) continue;
8469e24d8a0SVaclav Hapla     ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
8479e24d8a0SVaclav Hapla     ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
848ea211068SVaclav 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);}
849ea211068SVaclav 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])) {
850f80536cbSVaclav Hapla       PetscInt masterCone[2];
851f80536cbSVaclav Hapla       /* Translate these two cone points back to leave numbering */
852f80536cbSVaclav Hapla       for (c = 0; c < 2; c++) {
853def636c8SVaclav Hapla         if (leavesRanks[p][c] == rank) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "this should never happen - remote rank of point %D is the same rank",leavesRanks[p][c]);
854f80536cbSVaclav Hapla         /* Find index of rank leavesRanks[p][c] among remote ranks */
855f80536cbSVaclav Hapla         /* No need for PetscMPIIntCast because these integers were originally cast from PetscMPIInt. */
856f80536cbSVaclav Hapla         ierr = PetscFindMPIInt((PetscMPIInt)leavesRanks[p][c], nranks, ranks, &r);CHKERRQ(ierr);
857f80536cbSVaclav Hapla         if (PetscUnlikely(r < 0)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "this should never happen - rank %D not found among remote ranks",leavesRanks[p][c]);
858f80536cbSVaclav Hapla         /* Find point leaves[p][c] among remote points aimed at rank leavesRanks[p][c] */
859f80536cbSVaclav Hapla         o = roffset[r];
860f80536cbSVaclav Hapla         n = roffset[r+1] - o;
861f80536cbSVaclav Hapla         ierr = PetscFindInt(leaves[p][c], n, &rremote1[o], &ind0);CHKERRQ(ierr);
862adeface4SVaclav Hapla         if (PetscUnlikely(ind0 < 0)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "No cone point of %D is connected to (%D, %D) - it seems there is missing connection in point SF!",p,ranks[r],leaves[p][c]);
863f80536cbSVaclav Hapla         /* Get the corresponding local point */
864f80536cbSVaclav Hapla         masterCone[c] = rmine1[o+ind0];CHKERRQ(ierr);
865f80536cbSVaclav Hapla       }
866f80536cbSVaclav Hapla       if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d]  %D: masterCone=[%D %D]\n", rank, p, masterCone[0], masterCone[1]);CHKERRQ(ierr);}
867f80536cbSVaclav Hapla       /* Vaclav's note: Here we only compare first 2 points of the cone. Full cone size would lead to stronger self-checking. */
868f80536cbSVaclav Hapla       ierr = DMPlexOrientCell(dm, p, 2, masterCone);CHKERRQ(ierr);
86923aaf07bSVaclav Hapla     }
8702e745bdaSMatthew G. Knepley   }
871adeface4SVaclav Hapla #if defined(PETSC_USE_DEBUG)
872adeface4SVaclav Hapla   ierr = DMViewFromOptions(dm, NULL, "-after_fix_dm_view");CHKERRQ(ierr);
873adeface4SVaclav Hapla   for (r = 0; r < nleaves; ++r) {
874adeface4SVaclav Hapla     p = locals[r];
875adeface4SVaclav Hapla     ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
876adeface4SVaclav Hapla     if (!coneSize) continue;
877adeface4SVaclav Hapla     ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
878adeface4SVaclav Hapla     for (c = 0; c < 2; c++) {
879adeface4SVaclav Hapla       ierr = PetscFindInt(cone[c], nleaves, locals, &ind0);CHKERRQ(ierr);
880adeface4SVaclav Hapla       if (ind0 < 0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point SF contains %D but is missing its cone point cone[%D] = %D!", p, c, cone[c]);
881adeface4SVaclav Hapla       if (leaves[p][c] != remotes[ind0].index || leavesRanks[p][c] != remotes[ind0].rank) {
882adeface4SVaclav Hapla         if (leavesRanks[p][c] == rank) {
883adeface4SVaclav Hapla           PetscInt ind1;
884adeface4SVaclav Hapla           ierr = PetscFindInt(leaves[p][c], nleaves, locals, &ind1);CHKERRQ(ierr);
885adeface4SVaclav Hapla           if (ind1 < 0) {
886adeface4SVaclav Hapla             SETERRQ8(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D = locals[%d]: cone[%D]=%D --> (%D, %D) differs from the enforced (%D, %D). The latter was not even found among the local SF points - it is probably broken!", p, r, c, cone[c], remotes[ind0].rank, remotes[ind0].index, leavesRanks[p][c], leaves[p][c]);
887adeface4SVaclav Hapla           } else {
888adeface4SVaclav Hapla             SETERRQ9(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D = locals[%d]: cone[%D]=%D --> (%D, %D) differs from the enforced %D --> (%D, %D). Is the algorithm above or the point SF broken?", p, r, c, cone[c], remotes[ind0].rank, remotes[ind0].index, leaves[p][c], remotes[ind1].rank, remotes[ind1].index);
889adeface4SVaclav Hapla           }
890adeface4SVaclav Hapla         } else {
891adeface4SVaclav Hapla           SETERRQ8(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D = locals[%d]: cone[%D]=%D --> (%D, %D) differs from the enforced (%D, %D). Is the algorithm above or the point SF broken?", p, r, c, cone[c], remotes[ind0].rank, remotes[ind0].index, leavesRanks[p][c], leaves[p][c]);
892adeface4SVaclav Hapla         }
893adeface4SVaclav Hapla       }
894adeface4SVaclav Hapla     }
895adeface4SVaclav Hapla   }
896adeface4SVaclav Hapla #endif
8972e745bdaSMatthew G. Knepley   if (debug) {ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);}
8988a650c75SVaclav Hapla   ierr = PetscFree4(roots, leaves, rootsRanks, leavesRanks);CHKERRQ(ierr);
8997c7bb832SVaclav Hapla   ierr = PetscFree2(rmine1, rremote1);CHKERRQ(ierr);
9002e745bdaSMatthew G. Knepley   PetscFunctionReturn(0);
9012e745bdaSMatthew G. Knepley }
9022e745bdaSMatthew G. Knepley 
9032e72742cSMatthew G. Knepley static PetscErrorCode IntArrayViewFromOptions(MPI_Comm comm, const char opt[], const char name[], const char idxname[], const char valname[], PetscInt n, const PetscInt a[])
9047bffde78SMichael Lange {
9052e72742cSMatthew G. Knepley   PetscInt       idx;
9062e72742cSMatthew G. Knepley   PetscMPIInt    rank;
9072e72742cSMatthew G. Knepley   PetscBool      flg;
9087bffde78SMichael Lange   PetscErrorCode ierr;
9097bffde78SMichael Lange 
9107bffde78SMichael Lange   PetscFunctionBegin;
9112e72742cSMatthew G. Knepley   ierr = PetscOptionsHasName(NULL, NULL, opt, &flg);CHKERRQ(ierr);
9122e72742cSMatthew G. Knepley   if (!flg) PetscFunctionReturn(0);
9132e72742cSMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
9142e72742cSMatthew G. Knepley   ierr = PetscSynchronizedPrintf(comm, "[%d]%s:\n", rank, name);CHKERRQ(ierr);
9152e72742cSMatthew G. Knepley   for (idx = 0; idx < n; ++idx) {ierr = PetscSynchronizedPrintf(comm, "[%d]%s %D %s %D\n", rank, idxname, idx, valname, a[idx]);CHKERRQ(ierr);}
9162e72742cSMatthew G. Knepley   ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);
9172e72742cSMatthew G. Knepley   PetscFunctionReturn(0);
9182e72742cSMatthew G. Knepley }
9192e72742cSMatthew G. Knepley 
9202e72742cSMatthew G. Knepley static PetscErrorCode SFNodeArrayViewFromOptions(MPI_Comm comm, const char opt[], const char name[], const char idxname[], PetscInt n, const PetscSFNode a[])
9212e72742cSMatthew G. Knepley {
9222e72742cSMatthew G. Knepley   PetscInt       idx;
9232e72742cSMatthew G. Knepley   PetscMPIInt    rank;
9242e72742cSMatthew G. Knepley   PetscBool      flg;
9252e72742cSMatthew G. Knepley   PetscErrorCode ierr;
9262e72742cSMatthew G. Knepley 
9272e72742cSMatthew G. Knepley   PetscFunctionBegin;
9282e72742cSMatthew G. Knepley   ierr = PetscOptionsHasName(NULL, NULL, opt, &flg);CHKERRQ(ierr);
9292e72742cSMatthew G. Knepley   if (!flg) PetscFunctionReturn(0);
9302e72742cSMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
9312e72742cSMatthew G. Knepley   ierr = PetscSynchronizedPrintf(comm, "[%d]%s:\n", rank, name);CHKERRQ(ierr);
9322e72742cSMatthew G. Knepley   if (idxname) {
9332e72742cSMatthew G. Knepley     for (idx = 0; idx < n; ++idx) {ierr = PetscSynchronizedPrintf(comm, "[%d]%s %D rank %D index %D\n", rank, idxname, idx, a[idx].rank, a[idx].index);CHKERRQ(ierr);}
9342e72742cSMatthew G. Knepley   } else {
9352e72742cSMatthew G. Knepley     for (idx = 0; idx < n; ++idx) {ierr = PetscSynchronizedPrintf(comm, "[%d]rank %D index %D\n", rank, a[idx].rank, a[idx].index);CHKERRQ(ierr);}
9362e72742cSMatthew G. Knepley   }
9372e72742cSMatthew G. Knepley   ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);
9382e72742cSMatthew G. Knepley   PetscFunctionReturn(0);
9392e72742cSMatthew G. Knepley }
9402e72742cSMatthew G. Knepley 
9412e72742cSMatthew G. Knepley static PetscErrorCode DMPlexMapToLocalPoint(PetscHMapIJ roothash, const PetscInt localPoints[], PetscMPIInt rank, PetscSFNode remotePoint, PetscInt *localPoint)
9422e72742cSMatthew G. Knepley {
9432e72742cSMatthew G. Knepley   PetscErrorCode ierr;
9442e72742cSMatthew G. Knepley 
9452e72742cSMatthew G. Knepley   PetscFunctionBegin;
9462e72742cSMatthew G. Knepley   if (remotePoint.rank == rank) {
9472e72742cSMatthew G. Knepley     *localPoint = remotePoint.index;
9482e72742cSMatthew G. Knepley   } else {
9492e72742cSMatthew G. Knepley     PetscHashIJKey key;
9502e72742cSMatthew G. Knepley     PetscInt       root;
9512e72742cSMatthew G. Knepley 
9522e72742cSMatthew G. Knepley     key.i = remotePoint.index;
9532e72742cSMatthew G. Knepley     key.j = remotePoint.rank;
9542e72742cSMatthew G. Knepley     ierr = PetscHMapIJGet(roothash, key, &root);CHKERRQ(ierr);
9552e72742cSMatthew G. Knepley     if (root >= 0) {
9562e72742cSMatthew G. Knepley       *localPoint = localPoints[root];
9572e72742cSMatthew G. Knepley     } else PetscFunctionReturn(1);
9582e72742cSMatthew G. Knepley   }
9592e72742cSMatthew G. Knepley   PetscFunctionReturn(0);
9602e72742cSMatthew G. Knepley }
9612e72742cSMatthew G. Knepley 
9622e72742cSMatthew G. Knepley /*@
9632e72742cSMatthew G. Knepley   DMPlexInterpolatePointSF - Insert interpolated points in the overlap into the PointSF in parallel, following local interpolation
9642e72742cSMatthew G. Knepley 
9652e72742cSMatthew G. Knepley   Collective on DM
9662e72742cSMatthew G. Knepley 
9672e72742cSMatthew G. Knepley   Input Parameters:
9682e72742cSMatthew G. Knepley + dm      - The interpolated DM
9692e72742cSMatthew G. Knepley - pointSF - The initial SF without interpolated points
9702e72742cSMatthew G. Knepley 
9712e72742cSMatthew G. Knepley   Output Parameter:
9722e72742cSMatthew G. Knepley . pointSF - The SF including interpolated points
9732e72742cSMatthew G. Knepley 
9742e72742cSMatthew G. Knepley   Level: intermediate
9752e72742cSMatthew G. Knepley 
9762e72742cSMatthew G. Knepley    Note: All debugging for this process can be turned on with the options: -dm_interp_pre_view -petscsf_interp_pre_view -petscsection_interp_candidate_view -petscsection_interp_candidate_remote_view -petscsection_interp_claim_view -petscsf_interp_pre_view -dmplex_interp_debug
9772e72742cSMatthew G. Knepley 
9782e72742cSMatthew G. Knepley .keywords: mesh
9792e72742cSMatthew G. Knepley .seealso: DMPlexInterpolate(), DMPlexUninterpolate()
9802e72742cSMatthew G. Knepley @*/
981e53487d3SMatthew G. Knepley PetscErrorCode DMPlexInterpolatePointSF(DM dm, PetscSF pointSF)
9822e72742cSMatthew G. Knepley {
9832e72742cSMatthew G. Knepley   /*
9842e72742cSMatthew G. Knepley        Okay, the algorithm is:
9852e72742cSMatthew G. Knepley          - Take each point in the overlap (root)
9862e72742cSMatthew G. Knepley          - Look at the neighboring points in the overlap (candidates)
9872e72742cSMatthew G. Knepley          - Send these candidate points to neighbors
9882e72742cSMatthew G. Knepley          - Neighbor checks for edge between root and candidate
9892e72742cSMatthew G. Knepley          - If edge is found, it replaces candidate point with edge point
9902e72742cSMatthew G. Knepley          - Send back the overwritten candidates (claims)
9912e72742cSMatthew G. Knepley          - Original guy checks for edges, different from original candidate, and gets its own edge
9922e72742cSMatthew G. Knepley          - This pair is put into SF
9932e72742cSMatthew G. Knepley 
9942e72742cSMatthew G. Knepley        We need a new algorithm that tolerates groups larger than 2.
9952e72742cSMatthew G. Knepley          - Take each point in the overlap (root)
9962e72742cSMatthew G. Knepley          - Find all collections of points in the overlap which make faces (do early join)
9972e72742cSMatthew G. Knepley          - Send collections as candidates (add size as first number)
998*66aa2a39SMatthew G. Knepley            - Make sure to send collection to all owners of all overlap points in collection
9992e72742cSMatthew G. Knepley          - Neighbor check for face in collections
10002e72742cSMatthew G. Knepley          - If face is found, it replaces candidate point with face point
10012e72742cSMatthew G. Knepley          - Send back the overwritten candidates (claims)
10022e72742cSMatthew G. Knepley          - Original guy checks for faces, different from original candidate, and gets its own face
10032e72742cSMatthew G. Knepley          - This pair is put into SF
10042e72742cSMatthew G. Knepley   */
10052e72742cSMatthew G. Knepley   PetscHMapI         leafhash;
10062e72742cSMatthew G. Knepley   PetscHMapIJ        roothash;
10072e72742cSMatthew G. Knepley   const PetscInt    *localPoints, *rootdegree;
10082e72742cSMatthew G. Knepley   const PetscSFNode *remotePoints;
10092e72742cSMatthew G. Knepley   PetscSFNode       *candidates, *candidatesRemote, *claims;
10102e72742cSMatthew G. Knepley   PetscSection       candidateSection, candidateSectionRemote, claimSection;
10112e72742cSMatthew G. Knepley   PetscInt           numLeaves, l, numRoots, r, candidatesSize, candidatesRemoteSize;
10122e72742cSMatthew G. Knepley   PetscMPIInt        size, rank;
10132e72742cSMatthew G. Knepley   PetscHashIJKey     key;
10142e72742cSMatthew G. Knepley   PetscBool          debug = PETSC_FALSE;
10152e72742cSMatthew G. Knepley   PetscErrorCode     ierr;
10162e72742cSMatthew G. Knepley 
10172e72742cSMatthew G. Knepley   PetscFunctionBegin;
10182e72742cSMatthew G. Knepley   ierr = PetscOptionsHasName(NULL, ((PetscObject) dm)->prefix, "-dmplex_interp_debug", &debug);CHKERRQ(ierr);
10199852e123SBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr);
10207bffde78SMichael Lange   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
10217bffde78SMichael Lange   ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr);
10229852e123SBarry Smith   if (size < 2 || numRoots < 0) PetscFunctionReturn(0);
10232e72742cSMatthew G. Knepley   ierr = PetscObjectViewFromOptions((PetscObject) dm, NULL, "-dm_interp_pre_view");CHKERRQ(ierr);
10242e72742cSMatthew G. Knepley   ierr = PetscObjectViewFromOptions((PetscObject) pointSF, NULL, "-petscsf_interp_pre_view");CHKERRQ(ierr);
102525afeb17SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr);
10267bffde78SMichael Lange   /* Build hashes of points in the SF for efficient lookup */
1027e8f14785SLisandro Dalcin   ierr = PetscHMapICreate(&leafhash);CHKERRQ(ierr);
1028e8f14785SLisandro Dalcin   ierr = PetscHMapIJCreate(&roothash);CHKERRQ(ierr);
10292e72742cSMatthew G. Knepley   for (l = 0; l < numLeaves; ++l) {
10302e72742cSMatthew G. Knepley     key.i = remotePoints[l].index;
10312e72742cSMatthew G. Knepley     key.j = remotePoints[l].rank;
10322e72742cSMatthew G. Knepley     ierr = PetscHMapISet(leafhash, localPoints[l], l);CHKERRQ(ierr);
10332e72742cSMatthew G. Knepley     ierr = PetscHMapIJSet(roothash, key, l);CHKERRQ(ierr);
10347bffde78SMichael Lange   }
1035*66aa2a39SMatthew G. Knepley   /* Compute root degree to identify shared points */
10362e72742cSMatthew G. Knepley   ierr = PetscSFComputeDegreeBegin(pointSF, &rootdegree);CHKERRQ(ierr);
10372e72742cSMatthew G. Knepley   ierr = PetscSFComputeDegreeEnd(pointSF, &rootdegree);CHKERRQ(ierr);
10382e72742cSMatthew G. Knepley   ierr = IntArrayViewFromOptions(PetscObjectComm((PetscObject) dm), "-interp_root_degree_view", "Root degree", "point", "degree", numRoots, rootdegree);CHKERRQ(ierr);
1039*66aa2a39SMatthew G. Knepley   /* Build a section / SFNode array of candidate points (face bd points) in the cone(support(leaf)),
1040*66aa2a39SMatthew G. Knepley      where each candidate is defined by a set of remote points (roots) for the other points that define the face. */
10417bffde78SMichael Lange   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSection);CHKERRQ(ierr);
10427bffde78SMichael Lange   ierr = PetscSectionSetChart(candidateSection, 0, numRoots);CHKERRQ(ierr);
10437bffde78SMichael Lange   {
1044*66aa2a39SMatthew G. Knepley     PetscHMapIJ facehash;
10452e72742cSMatthew G. Knepley 
1046*66aa2a39SMatthew G. Knepley     ierr = PetscHMapIJCreate(&facehash);CHKERRQ(ierr);
10472e72742cSMatthew G. Knepley     for (l = 0; l < numLeaves; ++l) {
10482e72742cSMatthew G. Knepley       const PetscInt    localPoint = localPoints[l];
10492e72742cSMatthew G. Knepley       const PetscInt   *support;
10502e72742cSMatthew G. Knepley       PetscInt          supportSize, s;
10512e72742cSMatthew G. Knepley 
10522e72742cSMatthew G. Knepley       if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]  Checking local point %D\n", rank, localPoint);CHKERRQ(ierr);}
10532e72742cSMatthew G. Knepley       ierr = DMPlexGetSupportSize(dm, localPoint, &supportSize);CHKERRQ(ierr);
10542e72742cSMatthew G. Knepley       ierr = DMPlexGetSupport(dm, localPoint, &support);CHKERRQ(ierr);
10552e72742cSMatthew G. Knepley       for (s = 0; s < supportSize; ++s) {
10562e72742cSMatthew G. Knepley         const PetscInt  face = support[s];
10572e72742cSMatthew G. Knepley         const PetscInt *cone;
10582e72742cSMatthew G. Knepley         PetscInt        coneSize, c, f, root;
10592e72742cSMatthew G. Knepley         PetscBool       isFace = PETSC_TRUE;
10602e72742cSMatthew G. Knepley 
10612e72742cSMatthew G. Knepley         /* Only add face once */
10622e72742cSMatthew G. Knepley         if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]    Support point %D\n", rank, face);CHKERRQ(ierr);}
1063*66aa2a39SMatthew G. Knepley         key.i = localPoint;
1064*66aa2a39SMatthew G. Knepley         key.j = face;
1065*66aa2a39SMatthew G. Knepley         ierr = PetscHMapIJGet(facehash, key, &f);CHKERRQ(ierr);
10662e72742cSMatthew G. Knepley         if (f >= 0) continue;
10672e72742cSMatthew G. Knepley         ierr = DMPlexGetConeSize(dm, face, &coneSize);CHKERRQ(ierr);
10682e72742cSMatthew G. Knepley         ierr = DMPlexGetCone(dm, face, &cone);CHKERRQ(ierr);
10692e72742cSMatthew G. Knepley         /* If a cone point does not map to leaves on any proc, then do not put face in SF */
10702e72742cSMatthew G. Knepley         for (c = 0; c < coneSize; ++c) {
10712e72742cSMatthew G. Knepley           if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]      Cone point %D\n", rank, cone[c]);CHKERRQ(ierr);}
10722e72742cSMatthew G. Knepley           ierr = PetscHMapIGet(leafhash, cone[c], &root);CHKERRQ(ierr);
10732e72742cSMatthew G. Knepley           if (!rootdegree[cone[c]] && (root < 0)) {isFace = PETSC_FALSE; break;}
10742e72742cSMatthew G. Knepley         }
10752e72742cSMatthew G. Knepley         if (isFace) {
1076*66aa2a39SMatthew G. Knepley           if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]    Found shared face %D\n", rank, face);CHKERRQ(ierr);}
1077*66aa2a39SMatthew G. Knepley           ierr = PetscHMapIJSet(facehash, key, l);CHKERRQ(ierr);
10782e72742cSMatthew G. Knepley           ierr = PetscSectionAddDof(candidateSection, localPoint, coneSize);CHKERRQ(ierr);
10797bffde78SMichael Lange         }
10807bffde78SMichael Lange       }
10812e72742cSMatthew G. Knepley     }
10822e72742cSMatthew G. Knepley     if (debug) {ierr = PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);CHKERRQ(ierr);}
1083*66aa2a39SMatthew G. Knepley     ierr = PetscHMapIJClear(facehash);CHKERRQ(ierr);
10847bffde78SMichael Lange     ierr = PetscSectionSetUp(candidateSection);CHKERRQ(ierr);
10857bffde78SMichael Lange     ierr = PetscSectionGetStorageSize(candidateSection, &candidatesSize);CHKERRQ(ierr);
10867bffde78SMichael Lange     ierr = PetscMalloc1(candidatesSize, &candidates);CHKERRQ(ierr);
10872e72742cSMatthew G. Knepley     for (l = 0; l < numLeaves; ++l) {
10882e72742cSMatthew G. Knepley       const PetscInt    localPoint = localPoints[l];
10892e72742cSMatthew G. Knepley       const PetscInt   *support;
10902e72742cSMatthew G. Knepley       PetscInt          supportSize, s, offset, idx = 0;
10912e72742cSMatthew G. Knepley 
10922e72742cSMatthew G. Knepley       if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]  Checking local point %D\n", rank, localPoint);CHKERRQ(ierr);}
10932e72742cSMatthew G. Knepley       ierr = PetscSectionGetOffset(candidateSection, localPoint, &offset);CHKERRQ(ierr);
10942e72742cSMatthew G. Knepley       ierr = DMPlexGetSupportSize(dm, localPoint, &supportSize);CHKERRQ(ierr);
10952e72742cSMatthew G. Knepley       ierr = DMPlexGetSupport(dm, localPoint, &support);CHKERRQ(ierr);
10962e72742cSMatthew G. Knepley       for (s = 0; s < supportSize; ++s) {
10972e72742cSMatthew G. Knepley         const PetscInt  face = support[s];
10982e72742cSMatthew G. Knepley         const PetscInt *cone;
10992e72742cSMatthew G. Knepley         PetscInt        coneSize, c, f, root;
11002e72742cSMatthew G. Knepley         PetscBool       isFace = PETSC_TRUE;
11012e72742cSMatthew G. Knepley 
11022e72742cSMatthew G. Knepley         /* Only add face once */
11032e72742cSMatthew G. Knepley         if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]    Support point %D\n", rank, face);CHKERRQ(ierr);}
1104*66aa2a39SMatthew G. Knepley         key.i = localPoint;
1105*66aa2a39SMatthew G. Knepley         key.j = face;
1106*66aa2a39SMatthew G. Knepley         ierr = PetscHMapIJGet(facehash, key, &f);CHKERRQ(ierr);
11072e72742cSMatthew G. Knepley         if (f >= 0) continue;
11082e72742cSMatthew G. Knepley         ierr = DMPlexGetConeSize(dm, face, &coneSize);CHKERRQ(ierr);
11092e72742cSMatthew G. Knepley         ierr = DMPlexGetCone(dm, face, &cone);CHKERRQ(ierr);
11102e72742cSMatthew G. Knepley         /* If a cone point does not map to leaves on any proc, then do not put face in SF */
11112e72742cSMatthew G. Knepley         for (c = 0; c < coneSize; ++c) {
11122e72742cSMatthew G. Knepley           if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]      Cone point %D\n", rank, cone[c]);CHKERRQ(ierr);}
11132e72742cSMatthew G. Knepley           ierr = PetscHMapIGet(leafhash, cone[c], &root);CHKERRQ(ierr);
11142e72742cSMatthew G. Knepley           if (!rootdegree[cone[c]] && (root < 0)) {isFace = PETSC_FALSE; break;}
11152e72742cSMatthew G. Knepley         }
11162e72742cSMatthew G. Knepley         if (isFace) {
1117*66aa2a39SMatthew G. Knepley           if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]    Adding shared face %D at idx %D\n", rank, face, idx);CHKERRQ(ierr);}
1118*66aa2a39SMatthew G. Knepley           ierr = PetscHMapIJSet(facehash, key, l);CHKERRQ(ierr);
11192e72742cSMatthew G. Knepley           candidates[offset+idx].rank    = -1;
11202e72742cSMatthew G. Knepley           candidates[offset+idx++].index = coneSize-1;
11212e72742cSMatthew G. Knepley           for (c = 0; c < coneSize; ++c) {
11222e72742cSMatthew G. Knepley             if (cone[c] == localPoint) continue;
11232e72742cSMatthew G. Knepley             if (rootdegree[cone[c]]) {
11242e72742cSMatthew G. Knepley               candidates[offset+idx].rank    = rank;
11252e72742cSMatthew G. Knepley               candidates[offset+idx++].index = cone[c];
11262e72742cSMatthew G. Knepley             } else {
11272e72742cSMatthew G. Knepley               ierr = PetscHMapIGet(leafhash, cone[c], &root);CHKERRQ(ierr);
11282e72742cSMatthew G. Knepley               if (root < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot locate local point %D in SF", cone[c]);
11292e72742cSMatthew G. Knepley               candidates[offset+idx++] = remotePoints[root];
11307bffde78SMichael Lange             }
11317bffde78SMichael Lange           }
11322e72742cSMatthew G. Knepley         }
11332e72742cSMatthew G. Knepley       }
11342e72742cSMatthew G. Knepley     }
11352e72742cSMatthew G. Knepley     if (debug) {ierr = PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);CHKERRQ(ierr);}
1136*66aa2a39SMatthew G. Knepley     ierr = PetscHMapIJDestroy(&facehash);CHKERRQ(ierr);
11372e72742cSMatthew G. Knepley     ierr = PetscObjectViewFromOptions((PetscObject) candidateSection, NULL, "-petscsection_interp_candidate_view");CHKERRQ(ierr);
11382e72742cSMatthew G. Knepley     ierr = SFNodeArrayViewFromOptions(PetscObjectComm((PetscObject) dm), "-petscsection_interp_candidate_view", "Candidates", NULL, candidatesSize, candidates);CHKERRQ(ierr);
11397bffde78SMichael Lange   }
11407bffde78SMichael Lange   /* Gather candidate section / array pair into the root partition via inverse(multi(pointSF)). */
11412e72742cSMatthew G. Knepley   /*   Note that this section is indexed by offsets into leaves, not by point number */
11427bffde78SMichael Lange   {
11437bffde78SMichael Lange     PetscSF   sfMulti, sfInverse, sfCandidates;
11447bffde78SMichael Lange     PetscInt *remoteOffsets;
11452e72742cSMatthew G. Knepley 
11467bffde78SMichael Lange     ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr);
11477bffde78SMichael Lange     ierr = PetscSFCreateInverseSF(sfMulti, &sfInverse);CHKERRQ(ierr);
11487bffde78SMichael Lange     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSectionRemote);CHKERRQ(ierr);
11497bffde78SMichael Lange     ierr = PetscSFDistributeSection(sfInverse, candidateSection, &remoteOffsets, candidateSectionRemote);CHKERRQ(ierr);
11507bffde78SMichael Lange     ierr = PetscSFCreateSectionSF(sfInverse, candidateSection, remoteOffsets, candidateSectionRemote, &sfCandidates);CHKERRQ(ierr);
11517bffde78SMichael Lange     ierr = PetscSectionGetStorageSize(candidateSectionRemote, &candidatesRemoteSize);CHKERRQ(ierr);
11527bffde78SMichael Lange     ierr = PetscMalloc1(candidatesRemoteSize, &candidatesRemote);CHKERRQ(ierr);
11537bffde78SMichael Lange     ierr = PetscSFBcastBegin(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr);
11547bffde78SMichael Lange     ierr = PetscSFBcastEnd(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr);
11557bffde78SMichael Lange     ierr = PetscSFDestroy(&sfInverse);CHKERRQ(ierr);
11567bffde78SMichael Lange     ierr = PetscSFDestroy(&sfCandidates);CHKERRQ(ierr);
11577bffde78SMichael Lange     ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
11582e72742cSMatthew G. Knepley 
11592e72742cSMatthew G. Knepley     ierr = PetscObjectViewFromOptions((PetscObject) candidateSectionRemote, NULL, "-petscsection_interp_candidate_remote_view");CHKERRQ(ierr);
11602e72742cSMatthew G. Knepley     ierr = SFNodeArrayViewFromOptions(PetscObjectComm((PetscObject) dm), "-petscsection_interp_candidate_remote_view", "Remote Candidates", NULL, candidatesRemoteSize, candidatesRemote);CHKERRQ(ierr);
11617bffde78SMichael Lange   }
11622e72742cSMatthew G. Knepley   /* */
11637bffde78SMichael Lange   {
11642e72742cSMatthew G. Knepley     PetscInt idx;
11652e72742cSMatthew G. Knepley     /* There is a section point for every leaf attached to a given root point */
11662e72742cSMatthew G. Knepley     for (r = 0, idx = 0; r < numRoots; ++r) {
11672e72742cSMatthew G. Knepley       PetscInt deg;
11682e72742cSMatthew G. Knepley       for (deg = 0; deg < rootdegree[r]; ++deg, ++idx) {
11692e72742cSMatthew G. Knepley         PetscInt offset, dof, d;
11702e72742cSMatthew G. Knepley 
11717bffde78SMichael Lange         ierr = PetscSectionGetDof(candidateSectionRemote, idx, &dof);CHKERRQ(ierr);
11727bffde78SMichael Lange         ierr = PetscSectionGetOffset(candidateSectionRemote, idx, &offset);CHKERRQ(ierr);
11732e72742cSMatthew G. Knepley         for (d = 0; d < dof; ++d) {
11742e72742cSMatthew G. Knepley           const PetscInt  sizeInd   = offset+d;
11752e72742cSMatthew G. Knepley           const PetscInt  numPoints = candidatesRemote[sizeInd].index;
11762e72742cSMatthew G. Knepley           const PetscInt *join      = NULL;
11772e72742cSMatthew G. Knepley           PetscInt        points[1024], p, joinSize;
11782e72742cSMatthew G. Knepley 
11792e72742cSMatthew G. Knepley           points[0] = r;
11802e72742cSMatthew G. Knepley           for (p = 0; p < numPoints; ++p) {
11812e72742cSMatthew G. Knepley             ierr = DMPlexMapToLocalPoint(roothash, localPoints, rank, candidatesRemote[offset+(++d)], &points[p+1]);
11822e72742cSMatthew G. Knepley             if (ierr) {d += numPoints-1 - p; break;} /* We got a point not in our overlap */
11832e72742cSMatthew G. Knepley             if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]  Checking local candidate %D\n", rank, points[p+1]);CHKERRQ(ierr);}
11847bffde78SMichael Lange           }
11852e72742cSMatthew G. Knepley           if (ierr) continue;
11862e72742cSMatthew G. Knepley           ierr = DMPlexGetJoin(dm, numPoints+1, points, &joinSize, &join);CHKERRQ(ierr);
11877bffde78SMichael Lange           if (joinSize == 1) {
11882e72742cSMatthew G. Knepley             if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]    Adding face %D at idx %D\n", rank, join[0], sizeInd);CHKERRQ(ierr);}
11892e72742cSMatthew G. Knepley             candidatesRemote[sizeInd].rank  = rank;
11902e72742cSMatthew G. Knepley             candidatesRemote[sizeInd].index = join[0];
11917bffde78SMichael Lange           }
11922e72742cSMatthew G. Knepley           ierr = DMPlexRestoreJoin(dm, numPoints+1, points, &joinSize, &join);CHKERRQ(ierr);
11937bffde78SMichael Lange         }
11947bffde78SMichael Lange       }
11957bffde78SMichael Lange     }
11962e72742cSMatthew G. Knepley     if (debug) {ierr = PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);CHKERRQ(ierr);}
11977bffde78SMichael Lange   }
11987bffde78SMichael Lange   /* Push claims back to receiver via the MultiSF and derive new pointSF mapping on receiver */
11997bffde78SMichael Lange   {
12007bffde78SMichael Lange     PetscSF         sfMulti, sfClaims, sfPointNew;
12017bffde78SMichael Lange     PetscSFNode    *remotePointsNew;
12022e72742cSMatthew G. Knepley     PetscHMapI      claimshash;
12032e72742cSMatthew G. Knepley     PetscInt       *remoteOffsets, *localPointsNew;
12042e72742cSMatthew G. Knepley     PetscInt        claimsSize, pStart, pEnd, root, numLocalNew, p, d;
12052e72742cSMatthew G. Knepley 
12067bffde78SMichael Lange     ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr);
12077bffde78SMichael Lange     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &claimSection);CHKERRQ(ierr);
12087bffde78SMichael Lange     ierr = PetscSFDistributeSection(sfMulti, candidateSectionRemote, &remoteOffsets, claimSection);CHKERRQ(ierr);
12097bffde78SMichael Lange     ierr = PetscSFCreateSectionSF(sfMulti, candidateSectionRemote, remoteOffsets, claimSection, &sfClaims);CHKERRQ(ierr);
12102e72742cSMatthew G. Knepley     ierr = PetscSectionGetStorageSize(claimSection, &claimsSize);CHKERRQ(ierr);
12112e72742cSMatthew G. Knepley     ierr = PetscMalloc1(claimsSize, &claims);CHKERRQ(ierr);
12127bffde78SMichael Lange     ierr = PetscSFBcastBegin(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr);
12137bffde78SMichael Lange     ierr = PetscSFBcastEnd(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr);
12147bffde78SMichael Lange     ierr = PetscSFDestroy(&sfClaims);CHKERRQ(ierr);
12157bffde78SMichael Lange     ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
12162e72742cSMatthew G. Knepley     ierr = PetscObjectViewFromOptions((PetscObject) claimSection, NULL, "-petscsection_interp_claim_view");CHKERRQ(ierr);
12172e72742cSMatthew G. Knepley     ierr = SFNodeArrayViewFromOptions(PetscObjectComm((PetscObject) dm), "-petscsection_interp_claim_view", "Claims", NULL, claimsSize, claims);CHKERRQ(ierr);
12187bffde78SMichael Lange     /* Walk the original section of local supports and add an SF entry for each updated item */
1219e8f14785SLisandro Dalcin     ierr = PetscHMapICreate(&claimshash);CHKERRQ(ierr);
12207bffde78SMichael Lange     for (p = 0; p < numRoots; ++p) {
12212e72742cSMatthew G. Knepley       PetscInt dof, offset;
12222e72742cSMatthew G. Knepley 
12232e72742cSMatthew G. Knepley       if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]  Checking root for claims %D\n", rank, p);CHKERRQ(ierr);}
12247bffde78SMichael Lange       ierr = PetscSectionGetDof(candidateSection, p, &dof);CHKERRQ(ierr);
12257bffde78SMichael Lange       ierr = PetscSectionGetOffset(candidateSection, p, &offset);CHKERRQ(ierr);
12262e72742cSMatthew G. Knepley       for (d = 0; d < dof;) {
12272e72742cSMatthew G. Knepley         if (claims[offset+d].rank >= 0) {
12282e72742cSMatthew G. Knepley           const PetscInt  faceInd   = offset+d;
12292e72742cSMatthew G. Knepley           const PetscInt  numPoints = candidates[faceInd].index;
12302e72742cSMatthew G. Knepley           const PetscInt *join      = NULL;
12312e72742cSMatthew G. Knepley           PetscInt        joinSize, points[1024], c;
12322e72742cSMatthew G. Knepley 
12332e72742cSMatthew G. Knepley           if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]    Found claim for remote point (%D, %D)\n", rank, claims[faceInd].rank, claims[faceInd].index);CHKERRQ(ierr);}
12342e72742cSMatthew G. Knepley           points[0] = p;
12352e72742cSMatthew G. Knepley           if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]      point %D\n", rank, points[0]);CHKERRQ(ierr);}
12362e72742cSMatthew G. Knepley           for (c = 0, ++d; c < numPoints; ++c, ++d) {
1237e8f14785SLisandro Dalcin             key.i = candidates[offset+d].index;
1238e8f14785SLisandro Dalcin             key.j = candidates[offset+d].rank;
1239e8f14785SLisandro Dalcin             ierr = PetscHMapIJGet(roothash, key, &root);CHKERRQ(ierr);
12402e72742cSMatthew G. Knepley             points[c+1] = localPoints[root];
12412e72742cSMatthew G. Knepley             if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]      point %D\n", rank, points[c+1]);CHKERRQ(ierr);}
12422e72742cSMatthew G. Knepley           }
12432e72742cSMatthew G. Knepley           ierr = DMPlexGetJoin(dm, numPoints+1, points, &joinSize, &join);CHKERRQ(ierr);
12442e72742cSMatthew G. Knepley           if (joinSize == 1) {
12452e72742cSMatthew G. Knepley             if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]    Found local face %D\n", rank, join[0]);CHKERRQ(ierr);}
12462e72742cSMatthew G. Knepley             ierr = PetscHMapISet(claimshash, join[0], faceInd);CHKERRQ(ierr);
12472e72742cSMatthew G. Knepley           }
12482e72742cSMatthew G. Knepley           ierr = DMPlexRestoreJoin(dm, numPoints+1, points, &joinSize, &join);CHKERRQ(ierr);
12492e72742cSMatthew G. Knepley         } else d += claims[offset+d].index+1;
12507bffde78SMichael Lange       }
12517bffde78SMichael Lange     }
12522e72742cSMatthew G. Knepley     if (debug) {ierr = PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);CHKERRQ(ierr);}
12537bffde78SMichael Lange     /* Create new pointSF from hashed claims */
1254e8f14785SLisandro Dalcin     ierr = PetscHMapIGetSize(claimshash, &numLocalNew);CHKERRQ(ierr);
12557bffde78SMichael Lange     ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
12567bffde78SMichael Lange     ierr = PetscMalloc1(numLeaves + numLocalNew, &localPointsNew);CHKERRQ(ierr);
12577bffde78SMichael Lange     ierr = PetscMalloc1(numLeaves + numLocalNew, &remotePointsNew);CHKERRQ(ierr);
12587bffde78SMichael Lange     for (p = 0; p < numLeaves; ++p) {
12597bffde78SMichael Lange       localPointsNew[p] = localPoints[p];
12607bffde78SMichael Lange       remotePointsNew[p].index = remotePoints[p].index;
12617bffde78SMichael Lange       remotePointsNew[p].rank  = remotePoints[p].rank;
12627bffde78SMichael Lange     }
1263f3190fc2SToby Isaac     p = numLeaves;
1264e8f14785SLisandro Dalcin     ierr = PetscHMapIGetKeys(claimshash, &p, localPointsNew);CHKERRQ(ierr);
1265f3190fc2SToby Isaac     ierr = PetscSortInt(numLocalNew, &localPointsNew[numLeaves]);CHKERRQ(ierr);
12667bffde78SMichael Lange     for (p = numLeaves; p < numLeaves + numLocalNew; ++p) {
12672e72742cSMatthew G. Knepley       PetscInt offset;
1268e8f14785SLisandro Dalcin       ierr = PetscHMapIGet(claimshash, localPointsNew[p], &offset);CHKERRQ(ierr);
12697bffde78SMichael Lange       remotePointsNew[p] = claims[offset];
12707bffde78SMichael Lange     }
12717bffde78SMichael Lange     ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), &sfPointNew);CHKERRQ(ierr);
12727bffde78SMichael Lange     ierr = PetscSFSetGraph(sfPointNew, pEnd-pStart, numLeaves+numLocalNew, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr);
12737bffde78SMichael Lange     ierr = DMSetPointSF(dm, sfPointNew);CHKERRQ(ierr);
12747bffde78SMichael Lange     ierr = PetscSFDestroy(&sfPointNew);CHKERRQ(ierr);
1275e8f14785SLisandro Dalcin     ierr = PetscHMapIDestroy(&claimshash);CHKERRQ(ierr);
12767bffde78SMichael Lange   }
1277e8f14785SLisandro Dalcin   ierr = PetscHMapIDestroy(&leafhash);CHKERRQ(ierr);
1278e8f14785SLisandro Dalcin   ierr = PetscHMapIJDestroy(&roothash);CHKERRQ(ierr);
12797bffde78SMichael Lange   ierr = PetscSectionDestroy(&candidateSection);CHKERRQ(ierr);
12807bffde78SMichael Lange   ierr = PetscSectionDestroy(&candidateSectionRemote);CHKERRQ(ierr);
12817bffde78SMichael Lange   ierr = PetscSectionDestroy(&claimSection);CHKERRQ(ierr);
12827bffde78SMichael Lange   ierr = PetscFree(candidates);CHKERRQ(ierr);
12837bffde78SMichael Lange   ierr = PetscFree(candidatesRemote);CHKERRQ(ierr);
12847bffde78SMichael Lange   ierr = PetscFree(claims);CHKERRQ(ierr);
128525afeb17SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr);
12867bffde78SMichael Lange   PetscFunctionReturn(0);
12877bffde78SMichael Lange }
12887bffde78SMichael Lange 
12890c795ddcSMatthew G. Knepley /*@C
129080330477SMatthew G. Knepley   DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc.
129180330477SMatthew G. Knepley 
129280330477SMatthew G. Knepley   Collective on DM
129380330477SMatthew G. Knepley 
1294e56d480eSMatthew G. Knepley   Input Parameters:
1295e56d480eSMatthew G. Knepley + dm - The DMPlex object with only cells and vertices
129610a67516SMatthew G. Knepley - dmInt - The interpolated DM
129780330477SMatthew G. Knepley 
129880330477SMatthew G. Knepley   Output Parameter:
12994e3744c5SMatthew G. Knepley . dmInt - The complete DMPlex object
130080330477SMatthew G. Knepley 
130180330477SMatthew G. Knepley   Level: intermediate
130280330477SMatthew G. Knepley 
130395452b02SPatrick Sanan   Notes:
130495452b02SPatrick Sanan     It does not copy over the coordinates.
130543eeeb2dSStefano Zampini 
130680330477SMatthew G. Knepley .keywords: mesh
130743eeeb2dSStefano Zampini .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates()
130880330477SMatthew G. Knepley @*/
13099f074e33SMatthew G Knepley PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt)
13109f074e33SMatthew G Knepley {
13119a5b13a2SMatthew G Knepley   DM             idm, odm = dm;
13127bffde78SMichael Lange   PetscSF        sfPoint;
13132e1b13c2SMatthew G. Knepley   PetscInt       depth, dim, d;
131410a67516SMatthew G. Knepley   const char    *name;
1315b325530aSVaclav Hapla   PetscBool      flg=PETSC_TRUE;
13169f074e33SMatthew G Knepley   PetscErrorCode ierr;
13179f074e33SMatthew G Knepley 
13189f074e33SMatthew G Knepley   PetscFunctionBegin;
131910a67516SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
132010a67516SMatthew G. Knepley   PetscValidPointer(dmInt, 2);
1321a72f3261SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr);
13222e1b13c2SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1323c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1324b21b8912SMatthew G. Knepley   if ((depth == dim) || (dim <= 1)) {
132576b791aaSMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
132676b791aaSMatthew G. Knepley     idm  = dm;
1327b21b8912SMatthew G. Knepley   } else {
13289a5b13a2SMatthew G Knepley     for (d = 1; d < dim; ++d) {
13299a5b13a2SMatthew G Knepley       /* Create interpolated mesh */
133010a67516SMatthew G. Knepley       ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr);
13319a5b13a2SMatthew G Knepley       ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr);
1332c73cfb54SMatthew G. Knepley       ierr = DMSetDimension(idm, dim);CHKERRQ(ierr);
13337bffde78SMichael Lange       if (depth > 0) {
13347bffde78SMichael Lange         ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr);
13357bffde78SMichael Lange         ierr = DMGetPointSF(odm, &sfPoint);CHKERRQ(ierr);
13362e72742cSMatthew G. Knepley         ierr = DMPlexInterpolatePointSF(idm, sfPoint);CHKERRQ(ierr);
13377bffde78SMichael Lange       }
13389a5b13a2SMatthew G Knepley       if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);}
13399a5b13a2SMatthew G Knepley       odm = idm;
13409f074e33SMatthew G Knepley     }
134110a67516SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) dm,  &name);CHKERRQ(ierr);
134210a67516SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) idm,  name);CHKERRQ(ierr);
134310a67516SMatthew G. Knepley     ierr = DMPlexCopyCoordinates(dm, idm);CHKERRQ(ierr);
134410a67516SMatthew G. Knepley     ierr = DMCopyLabels(dm, idm);CHKERRQ(ierr);
1345b325530aSVaclav Hapla     ierr = PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_interpolate_orient_interfaces", &flg, NULL);CHKERRQ(ierr);
1346b325530aSVaclav Hapla     if (flg) {ierr = DMPlexOrientInterface(idm);CHKERRQ(ierr);}
134784699f70SSatish Balay   }
134843eeeb2dSStefano Zampini   {
134943eeeb2dSStefano Zampini     PetscBool            isper;
135043eeeb2dSStefano Zampini     const PetscReal      *maxCell, *L;
135143eeeb2dSStefano Zampini     const DMBoundaryType *bd;
135243eeeb2dSStefano Zampini 
135343eeeb2dSStefano Zampini     ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr);
135443eeeb2dSStefano Zampini     ierr = DMSetPeriodicity(idm,isper,maxCell,L,bd);CHKERRQ(ierr);
135543eeeb2dSStefano Zampini   }
13569a5b13a2SMatthew G Knepley   *dmInt = idm;
1357a72f3261SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr);
13589f074e33SMatthew G Knepley   PetscFunctionReturn(0);
13599f074e33SMatthew G Knepley }
136007b0a1fcSMatthew G Knepley 
136180330477SMatthew G. Knepley /*@
136280330477SMatthew G. Knepley   DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices
136380330477SMatthew G. Knepley 
136480330477SMatthew G. Knepley   Collective on DM
136580330477SMatthew G. Knepley 
136680330477SMatthew G. Knepley   Input Parameter:
136780330477SMatthew G. Knepley . dmA - The DMPlex object with initial coordinates
136880330477SMatthew G. Knepley 
136980330477SMatthew G. Knepley   Output Parameter:
137080330477SMatthew G. Knepley . dmB - The DMPlex object with copied coordinates
137180330477SMatthew G. Knepley 
137280330477SMatthew G. Knepley   Level: intermediate
137380330477SMatthew G. Knepley 
137480330477SMatthew G. Knepley   Note: This is typically used when adding pieces other than vertices to a mesh
137580330477SMatthew G. Knepley 
137680330477SMatthew G. Knepley .keywords: mesh
137765f90189SJed Brown .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMGetCoordinateSection()
137880330477SMatthew G. Knepley @*/
137907b0a1fcSMatthew G Knepley PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB)
138007b0a1fcSMatthew G Knepley {
138107b0a1fcSMatthew G Knepley   Vec            coordinatesA, coordinatesB;
138243eeeb2dSStefano Zampini   VecType        vtype;
138307b0a1fcSMatthew G Knepley   PetscSection   coordSectionA, coordSectionB;
138407b0a1fcSMatthew G Knepley   PetscScalar   *coordsA, *coordsB;
13850bedd6beSMatthew G. Knepley   PetscInt       spaceDim, Nf, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d;
138643eeeb2dSStefano Zampini   PetscInt       cStartA, cEndA, cStartB, cEndB, cS, cE;
138743eeeb2dSStefano Zampini   PetscBool      lc = PETSC_FALSE;
138807b0a1fcSMatthew G Knepley   PetscErrorCode ierr;
138907b0a1fcSMatthew G Knepley 
139007b0a1fcSMatthew G Knepley   PetscFunctionBegin;
139143eeeb2dSStefano Zampini   PetscValidHeaderSpecific(dmA, DM_CLASSID, 1);
139243eeeb2dSStefano Zampini   PetscValidHeaderSpecific(dmB, DM_CLASSID, 2);
139376b791aaSMatthew G. Knepley   if (dmA == dmB) PetscFunctionReturn(0);
139407b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr);
139507b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr);
139607b0a1fcSMatthew 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);
139743eeeb2dSStefano Zampini   ierr = DMPlexGetHeightStratum(dmA, 0, &cStartA, &cEndA);CHKERRQ(ierr);
139843eeeb2dSStefano Zampini   ierr = DMPlexGetHeightStratum(dmB, 0, &cStartB, &cEndB);CHKERRQ(ierr);
139969d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr);
140069d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr);
1401972bc18aSToby Isaac   if (coordSectionA == coordSectionB) PetscFunctionReturn(0);
14020bedd6beSMatthew G. Knepley   ierr = PetscSectionGetNumFields(coordSectionA, &Nf);CHKERRQ(ierr);
14030bedd6beSMatthew G. Knepley   if (!Nf) PetscFunctionReturn(0);
14040bedd6beSMatthew G. Knepley   if (Nf > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of coordinate fields must be 1, not %D", Nf);
1405df26b574SMatthew G. Knepley   if (!coordSectionB) {
1406df26b574SMatthew G. Knepley     PetscInt dim;
1407df26b574SMatthew G. Knepley 
1408df26b574SMatthew G. Knepley     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) coordSectionA), &coordSectionB);CHKERRQ(ierr);
1409df26b574SMatthew G. Knepley     ierr = DMGetCoordinateDim(dmA, &dim);CHKERRQ(ierr);
1410df26b574SMatthew G. Knepley     ierr = DMSetCoordinateSection(dmB, dim, coordSectionB);CHKERRQ(ierr);
1411df26b574SMatthew G. Knepley     ierr = PetscObjectDereference((PetscObject) coordSectionB);CHKERRQ(ierr);
1412df26b574SMatthew G. Knepley   }
141307b0a1fcSMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr);
141407b0a1fcSMatthew G Knepley   ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr);
141507b0a1fcSMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr);
141643eeeb2dSStefano Zampini   ierr = PetscSectionGetChart(coordSectionA, &cS, &cE);CHKERRQ(ierr);
141743eeeb2dSStefano Zampini   if (cStartA <= cS && cS < cEndA) { /* localized coordinates */
1418367003a6SStefano Zampini     if ((cEndA-cStartA) != (cEndB-cStartB)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of cells in first DM %D != %D in the second DM", cEndA-cStartA, cEndB-cStartB);
141943eeeb2dSStefano Zampini     cS = cS - cStartA + cStartB;
142043eeeb2dSStefano Zampini     cE = vEndB;
142143eeeb2dSStefano Zampini     lc = PETSC_TRUE;
142243eeeb2dSStefano Zampini   } else {
142343eeeb2dSStefano Zampini     cS = vStartB;
142443eeeb2dSStefano Zampini     cE = vEndB;
142543eeeb2dSStefano Zampini   }
142643eeeb2dSStefano Zampini   ierr = PetscSectionSetChart(coordSectionB, cS, cE);CHKERRQ(ierr);
142707b0a1fcSMatthew G Knepley   for (v = vStartB; v < vEndB; ++v) {
142807b0a1fcSMatthew G Knepley     ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr);
142907b0a1fcSMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr);
143007b0a1fcSMatthew G Knepley   }
143143eeeb2dSStefano Zampini   if (lc) { /* localized coordinates */
143243eeeb2dSStefano Zampini     PetscInt c;
143343eeeb2dSStefano Zampini 
143443eeeb2dSStefano Zampini     for (c = cS-cStartB; c < cEndB-cStartB; c++) {
143543eeeb2dSStefano Zampini       PetscInt dof;
143643eeeb2dSStefano Zampini 
143743eeeb2dSStefano Zampini       ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr);
143843eeeb2dSStefano Zampini       ierr = PetscSectionSetDof(coordSectionB, c + cStartB, dof);CHKERRQ(ierr);
143943eeeb2dSStefano Zampini       ierr = PetscSectionSetFieldDof(coordSectionB, c + cStartB, 0, dof);CHKERRQ(ierr);
144043eeeb2dSStefano Zampini     }
144143eeeb2dSStefano Zampini   }
144207b0a1fcSMatthew G Knepley   ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr);
144307b0a1fcSMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr);
144407b0a1fcSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr);
14458b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinatesB);CHKERRQ(ierr);
144607b0a1fcSMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr);
144707b0a1fcSMatthew G Knepley   ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr);
144843eeeb2dSStefano Zampini   ierr = VecGetBlockSize(coordinatesA, &d);CHKERRQ(ierr);
144943eeeb2dSStefano Zampini   ierr = VecSetBlockSize(coordinatesB, d);CHKERRQ(ierr);
145043eeeb2dSStefano Zampini   ierr = VecGetType(coordinatesA, &vtype);CHKERRQ(ierr);
145143eeeb2dSStefano Zampini   ierr = VecSetType(coordinatesB, vtype);CHKERRQ(ierr);
145207b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr);
145307b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr);
145407b0a1fcSMatthew G Knepley   for (v = 0; v < vEndB-vStartB; ++v) {
145543eeeb2dSStefano Zampini     PetscInt offA, offB;
145643eeeb2dSStefano Zampini 
145743eeeb2dSStefano Zampini     ierr = PetscSectionGetOffset(coordSectionA, v + vStartA, &offA);CHKERRQ(ierr);
145843eeeb2dSStefano Zampini     ierr = PetscSectionGetOffset(coordSectionB, v + vStartB, &offB);CHKERRQ(ierr);
145907b0a1fcSMatthew G Knepley     for (d = 0; d < spaceDim; ++d) {
146043eeeb2dSStefano Zampini       coordsB[offB+d] = coordsA[offA+d];
146143eeeb2dSStefano Zampini     }
146243eeeb2dSStefano Zampini   }
146343eeeb2dSStefano Zampini   if (lc) { /* localized coordinates */
146443eeeb2dSStefano Zampini     PetscInt c;
146543eeeb2dSStefano Zampini 
146643eeeb2dSStefano Zampini     for (c = cS-cStartB; c < cEndB-cStartB; c++) {
146743eeeb2dSStefano Zampini       PetscInt dof, offA, offB;
146843eeeb2dSStefano Zampini 
146943eeeb2dSStefano Zampini       ierr = PetscSectionGetOffset(coordSectionA, c + cStartA, &offA);CHKERRQ(ierr);
147043eeeb2dSStefano Zampini       ierr = PetscSectionGetOffset(coordSectionB, c + cStartB, &offB);CHKERRQ(ierr);
147143eeeb2dSStefano Zampini       ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr);
147243eeeb2dSStefano Zampini       ierr = PetscMemcpy(coordsB + offB,coordsA + offA,dof*sizeof(*coordsB));CHKERRQ(ierr);
147307b0a1fcSMatthew G Knepley     }
147407b0a1fcSMatthew G Knepley   }
147507b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr);
147607b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr);
147707b0a1fcSMatthew G Knepley   ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr);
147807b0a1fcSMatthew G Knepley   ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr);
147907b0a1fcSMatthew G Knepley   PetscFunctionReturn(0);
148007b0a1fcSMatthew G Knepley }
14815c386225SMatthew G. Knepley 
14824e3744c5SMatthew G. Knepley /*@
14834e3744c5SMatthew G. Knepley   DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh
14844e3744c5SMatthew G. Knepley 
14854e3744c5SMatthew G. Knepley   Collective on DM
14864e3744c5SMatthew G. Knepley 
14874e3744c5SMatthew G. Knepley   Input Parameter:
14884e3744c5SMatthew G. Knepley . dm - The complete DMPlex object
14894e3744c5SMatthew G. Knepley 
14904e3744c5SMatthew G. Knepley   Output Parameter:
14914e3744c5SMatthew G. Knepley . dmUnint - The DMPlex object with only cells and vertices
14924e3744c5SMatthew G. Knepley 
14934e3744c5SMatthew G. Knepley   Level: intermediate
14944e3744c5SMatthew G. Knepley 
149595452b02SPatrick Sanan   Notes:
149695452b02SPatrick Sanan     It does not copy over the coordinates.
149743eeeb2dSStefano Zampini 
14984e3744c5SMatthew G. Knepley .keywords: mesh
149943eeeb2dSStefano Zampini .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates()
15004e3744c5SMatthew G. Knepley @*/
15014e3744c5SMatthew G. Knepley PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint)
15024e3744c5SMatthew G. Knepley {
15034e3744c5SMatthew G. Knepley   DM             udm;
1504c9f63434SStefano Zampini   PetscInt       dim, vStart, vEnd, cStart, cEnd, cMax, c, maxConeSize = 0, *cone;
15054e3744c5SMatthew G. Knepley   PetscErrorCode ierr;
15064e3744c5SMatthew G. Knepley 
15074e3744c5SMatthew G. Knepley   PetscFunctionBegin;
150843eeeb2dSStefano Zampini   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
150943eeeb2dSStefano Zampini   PetscValidPointer(dmUnint, 2);
1510c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
15114e3744c5SMatthew G. Knepley   if (dim <= 1) {
15124e3744c5SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
1513595d4abbSMatthew G. Knepley     *dmUnint = dm;
1514595d4abbSMatthew G. Knepley     PetscFunctionReturn(0);
15154e3744c5SMatthew G. Knepley   }
1516595d4abbSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1517595d4abbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1518c9f63434SStefano Zampini   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
15194e3744c5SMatthew G. Knepley   ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr);
15204e3744c5SMatthew G. Knepley   ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr);
1521c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(udm, dim);CHKERRQ(ierr);
15224e3744c5SMatthew G. Knepley   ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr);
15234e3744c5SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1524595d4abbSMatthew G. Knepley     PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
15254e3744c5SMatthew G. Knepley 
15264e3744c5SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
15274e3744c5SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
15284e3744c5SMatthew G. Knepley       const PetscInt p = closure[cl];
15294e3744c5SMatthew G. Knepley 
15304e3744c5SMatthew G. Knepley       if ((p >= vStart) && (p < vEnd)) ++coneSize;
15314e3744c5SMatthew G. Knepley     }
15324e3744c5SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
15334e3744c5SMatthew G. Knepley     ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr);
1534595d4abbSMatthew G. Knepley     maxConeSize = PetscMax(maxConeSize, coneSize);
15354e3744c5SMatthew G. Knepley   }
15364e3744c5SMatthew G. Knepley   ierr = DMSetUp(udm);CHKERRQ(ierr);
1537785e854fSJed Brown   ierr = PetscMalloc1(maxConeSize, &cone);CHKERRQ(ierr);
15384e3744c5SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1539595d4abbSMatthew G. Knepley     PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
15404e3744c5SMatthew G. Knepley 
15414e3744c5SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
15424e3744c5SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
15434e3744c5SMatthew G. Knepley       const PetscInt p = closure[cl];
15444e3744c5SMatthew G. Knepley 
15454e3744c5SMatthew G. Knepley       if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p;
15464e3744c5SMatthew G. Knepley     }
15474e3744c5SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
15484e3744c5SMatthew G. Knepley     ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr);
15494e3744c5SMatthew G. Knepley   }
15504e3744c5SMatthew G. Knepley   ierr = PetscFree(cone);CHKERRQ(ierr);
1551c9f63434SStefano Zampini   ierr = DMPlexSetHybridBounds(udm, cMax, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
15524e3744c5SMatthew G. Knepley   ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr);
15534e3744c5SMatthew G. Knepley   ierr = DMPlexStratify(udm);CHKERRQ(ierr);
15545c7de58aSMatthew G. Knepley   /* Reduce SF */
15555c7de58aSMatthew G. Knepley   {
15565c7de58aSMatthew G. Knepley     PetscSF            sfPoint, sfPointUn;
15575c7de58aSMatthew G. Knepley     const PetscSFNode *remotePoints;
15585c7de58aSMatthew G. Knepley     const PetscInt    *localPoints;
15595c7de58aSMatthew G. Knepley     PetscSFNode       *remotePointsUn;
15605c7de58aSMatthew G. Knepley     PetscInt          *localPointsUn;
15615c7de58aSMatthew G. Knepley     PetscInt           vEnd, numRoots, numLeaves, l;
15625c7de58aSMatthew G. Knepley     PetscInt           numLeavesUn = 0, n = 0;
15635c7de58aSMatthew G. Knepley     PetscErrorCode     ierr;
15645c7de58aSMatthew G. Knepley 
15655c7de58aSMatthew G. Knepley     /* Get original SF information */
15665c7de58aSMatthew G. Knepley     ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
15675c7de58aSMatthew G. Knepley     ierr = DMGetPointSF(udm, &sfPointUn);CHKERRQ(ierr);
15685c7de58aSMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, 0, NULL, &vEnd);CHKERRQ(ierr);
15695c7de58aSMatthew G. Knepley     ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr);
15705c7de58aSMatthew G. Knepley     /* Allocate space for cells and vertices */
15715c7de58aSMatthew G. Knepley     for (l = 0; l < numLeaves; ++l) if (localPoints[l] < vEnd) numLeavesUn++;
15725c7de58aSMatthew G. Knepley     /* Fill in leaves */
15735c7de58aSMatthew G. Knepley     if (vEnd >= 0) {
15745c7de58aSMatthew G. Knepley       ierr = PetscMalloc1(numLeavesUn, &remotePointsUn);CHKERRQ(ierr);
15755c7de58aSMatthew G. Knepley       ierr = PetscMalloc1(numLeavesUn, &localPointsUn);CHKERRQ(ierr);
15765c7de58aSMatthew G. Knepley       for (l = 0; l < numLeaves; l++) {
15775c7de58aSMatthew G. Knepley         if (localPoints[l] < vEnd) {
15785c7de58aSMatthew G. Knepley           localPointsUn[n]        = localPoints[l];
15795c7de58aSMatthew G. Knepley           remotePointsUn[n].rank  = remotePoints[l].rank;
15805c7de58aSMatthew G. Knepley           remotePointsUn[n].index = remotePoints[l].index;
15815c7de58aSMatthew G. Knepley           ++n;
15825c7de58aSMatthew G. Knepley         }
15835c7de58aSMatthew G. Knepley       }
15845c7de58aSMatthew G. Knepley       if (n != numLeavesUn) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent number of leaves %d != %d", n, numLeavesUn);
15855c7de58aSMatthew G. Knepley       ierr = PetscSFSetGraph(sfPointUn, vEnd, numLeavesUn, localPointsUn, PETSC_OWN_POINTER, remotePointsUn, PETSC_OWN_POINTER);CHKERRQ(ierr);
15865c7de58aSMatthew G. Knepley     }
15875c7de58aSMatthew G. Knepley   }
158843eeeb2dSStefano Zampini   {
158943eeeb2dSStefano Zampini     PetscBool            isper;
159043eeeb2dSStefano Zampini     const PetscReal      *maxCell, *L;
159143eeeb2dSStefano Zampini     const DMBoundaryType *bd;
159243eeeb2dSStefano Zampini 
159343eeeb2dSStefano Zampini     ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr);
159443eeeb2dSStefano Zampini     ierr = DMSetPeriodicity(udm,isper,maxCell,L,bd);CHKERRQ(ierr);
159543eeeb2dSStefano Zampini   }
159643eeeb2dSStefano Zampini 
15974e3744c5SMatthew G. Knepley   *dmUnint = udm;
15984e3744c5SMatthew G. Knepley   PetscFunctionReturn(0);
15994e3744c5SMatthew G. Knepley }
1600